ATLAS Offline Software
Loading...
Searching...
No Matches
LArBarrelPresamplerCalculator Class Reference

#include <LArBarrelPresamplerCalculator.h>

Inheritance diagram for LArBarrelPresamplerCalculator:
Collaboration diagram for LArBarrelPresamplerCalculator:

Public Member Functions

 LArBarrelPresamplerCalculator (const std::string &name, ISvcLocator *pSvcLocator)
 LArBarrelPresamplerCalculator (const LArBarrelPresamplerCalculator &)=delete
LArBarrelPresamplerCalculatoroperator= (const LArBarrelPresamplerCalculator &)=delete
virtual StatusCode initialize () override final
virtual StatusCode finalize () override final
virtual G4float OOTcut () const override final
virtual G4bool Process (const G4Step *a_step, std::vector< LArHitData > &hdata) const override final
virtual G4bool isInTime (G4double hitTime) const override final
virtual void initializeForSDCreation () override

Protected Attributes

Gaudi::Property< bool > m_BirksLaw {this, "BirksLaw", true}
Gaudi::Property< double > m_Birksk {this, "Birksk", 0.05832}
Gaudi::Property< double > m_OOTcut {this, "OOTcut", 300*CLHEP::ns}

Private Attributes

ServiceHandle< ILArBarrelPresamplerGeometrym_geometry {this, "GeometryCalculator", "LArBarrelPresamplerGeometry"}
Gaudi::Property< bool > m_IflCur {this, "IflCur", true}
Gaudi::Property< std::string > m_detectorName {this, "DetectorName", "LArMgr"}
Gaudi::Property< bool > m_testbeam {this, "isTestbeam", false}
const PsMapm_psmap {nullptr}
const LArG4BirksLawm_birksLaw {nullptr}
G4String m_volname {"LArMgr::LAr::Barrel::Presampler"}

Detailed Description

Definition at line 25 of file LArBarrelPresamplerCalculator.h.

Constructor & Destructor Documentation

◆ LArBarrelPresamplerCalculator() [1/2]

LArBarrelPresamplerCalculator::LArBarrelPresamplerCalculator ( const std::string & name,
ISvcLocator * pSvcLocator )

Definition at line 41 of file LArBarrelPresamplerCalculator.cxx.

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

◆ LArBarrelPresamplerCalculator() [2/2]

LArBarrelPresamplerCalculator::LArBarrelPresamplerCalculator ( const LArBarrelPresamplerCalculator & )
delete

Member Function Documentation

◆ finalize()

StatusCode LArBarrelPresamplerCalculator::finalize ( )
finaloverridevirtual

Definition at line 86 of file LArBarrelPresamplerCalculator.cxx.

87{
88 if (m_BirksLaw) delete m_birksLaw;
89 return StatusCode::SUCCESS;
90}
Gaudi::Property< bool > m_BirksLaw

◆ initialize()

StatusCode LArBarrelPresamplerCalculator::initialize ( )
finaloverridevirtual

Definition at line 46 of file LArBarrelPresamplerCalculator.cxx.

47{
48 ATH_MSG_DEBUG("LArBarrelPresamplerCalculator: Beginning initialization ");
49 // Initialize private members.
50
51 // Access source of detector parameters.
53
54 // Initialize the geometry calculator
55 ATH_CHECK(m_geometry.retrieve());
56
57 if (m_IflCur)
58 {
59 ATH_MSG_DEBUG(" LArBarrelPresamplerCalculator: start reading current maps");
61 }
62
63 // Get the out-of-time cut from the detector parameters routine.
64 m_OOTcut = parameters->GetValue("LArExpHallOutOfTimeCut"); //FIXME should be done by configurables
65
66 if (m_BirksLaw)
67 {
68 const double Birks_LAr_density = 1.396;
69 m_birksLaw = new LArG4BirksLaw(Birks_LAr_density,m_Birksk);
70 ATH_MSG_INFO("Use Birks_LAr_density="<<Birks_LAr_density<<", Birksk="<<(double)m_Birksk);
71 }
72 else
73 {
74 ATH_MSG_DEBUG(" LArBarrelPresamplerCalculator: Birks' law OFF");
75 }
76
77 if(m_detectorName.empty()) m_volname="LAr::Barrel::Presampler";
78 else m_volname=m_detectorName+"::LAr::Barrel::Presampler";
79
80 ATH_MSG_DEBUG("End of LArBarrelPresamplerCalculator initialization ");
81
82 return StatusCode::SUCCESS;
83}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_INFO(x)
#define ATH_MSG_DEBUG(x)
LArGeo::VDetectorParameters LArVG4DetectorParameters
Gaudi::Property< std::string > m_detectorName
ServiceHandle< ILArBarrelPresamplerGeometry > m_geometry
Gaudi::Property< double > m_OOTcut
Gaudi::Property< double > m_Birksk
static const VDetectorParameters * GetInstance()
static const PsMap * GetPsMap()
Definition PsMap.cxx:40

◆ initializeForSDCreation()

virtual void LArCalculatorSvcImp::initializeForSDCreation ( )
inlineoverridevirtualinherited

Reimplemented in LArBarrelCalculator, and LArHECWheelCalculator.

Definition at line 19 of file LArCalculatorSvcImp.h.

19{};

◆ isInTime()

virtual G4bool LArBarrelPresamplerCalculator::isInTime ( G4double hitTime) const
inlinefinaloverridevirtual

Definition at line 43 of file LArBarrelPresamplerCalculator.h.

44 {
45 return !(std::fabs(hitTime) > m_OOTcut);
46 }
float hitTime(const AFP_SIDSimHit &hit)

◆ OOTcut()

virtual G4float LArBarrelPresamplerCalculator::OOTcut ( ) const
inlinefinaloverridevirtual

Definition at line 38 of file LArBarrelPresamplerCalculator.h.

38{ return m_OOTcut; }

◆ operator=()

LArBarrelPresamplerCalculator & LArBarrelPresamplerCalculator::operator= ( const LArBarrelPresamplerCalculator & )
delete

◆ Process()

G4bool LArBarrelPresamplerCalculator::Process ( const G4Step * a_step,
std::vector< LArHitData > & hdata ) const
finaloverridevirtual

Definition at line 93 of file LArBarrelPresamplerCalculator.cxx.

94{
95 hdata.clear();
96 LArG4Identifier identifier2;
97
98 // check the Step content is non trivial
99 G4double thisStepEnergyDeposit = a_step->GetTotalEnergyDeposit() * a_step->GetTrack()->GetWeight();
100 G4double thisStepLength = a_step->GetStepLength() / Units::mm;
101 G4double dstep = .1*Units::mm; // length of punctual charge for Current Option
102
103#ifdef DEBUGSTEP
104 ATH_MSG_DEBUG( "****** LArBarrelPresamplerCalculator: Step energy,length "
105 << thisStepEnergyDeposit << " " << thisStepLength);
106#endif
107 if(thisStepEnergyDeposit <= 0. || thisStepLength <= 0.)
108 {
109 return false;
110 }
111
112 // Get Step geometrical parameters (first and end)
113 G4StepPoint *thisStepPoint = a_step->GetPreStepPoint();
114 G4StepPoint *thisStepBackPoint = a_step->GetPostStepPoint();
115 G4ThreeVector startPoint = thisStepPoint->GetPosition();
116 G4ThreeVector endPoint = thisStepBackPoint->GetPosition();
117
118#ifdef DEBUGSTEP
119 ATH_MSG_DEBUG(" Global beginning step position "
120 << startPoint.x() << " " << startPoint.y() << " "
121 << startPoint.z());
122 ATH_MSG_DEBUG(" Global end step position "
123 << endPoint.x() << " " << endPoint.y() << " "
124 << endPoint.z());
125#endif
126
127
128 // find transformation to go inside local half presampler tube volume
129
130 const G4NavigationHistory* g4navigation = thisStepPoint->GetTouchable()->GetHistory();
131 G4int ndep = g4navigation->GetDepth();
132 G4int idep = -999;
133
134#ifdef DEBUGSTEP
135 ATH_MSG_DEBUG(" Detector Name " << m_detectorName);
136#endif
137
138 const G4String presamplerName((m_detectorName.empty()) ? "LAr::Barrel::Presampler" :
139 m_detectorName+"::LAr::Barrel::Presampler");
140 for (G4int ii=0;ii<=ndep;ii++) {
141 G4VPhysicalVolume* v1 = g4navigation->GetVolume(ii);
142#ifdef DEBUGSTEP
143 ATH_MSG_DEBUG(" Level,VolumeName " << ii << " " << v1->GetName());
144#endif
145 // FIXME More efficient to find the comparison volume once and compare pointers?
146 if (v1->GetName()==presamplerName) idep=ii;
147 }
148
149#ifdef DEBUGSTEP
150 ATH_MSG_DEBUG(" idep = " << idep);
151#endif
152
153 if (idep<0) {
154 ATH_MSG_WARNING(" LArBarrelPresamplerCalculator::Process Presampler volume not found !!");
155 return false;
156 }
157
158 // transformation to go from global frame to LAr::Barrel::Presampler frame
159 const G4AffineTransform transformation = g4navigation->GetTransform(idep);
160
161 // step beginning and end in local frame
162 G4ThreeVector startPointinLocal =
163 transformation.TransformPoint(startPoint);
164 G4ThreeVector endPointinLocal =
165 transformation.TransformPoint (endPoint);
166
167#ifdef DEBUGSTEP
168 ATH_MSG_DEBUG(" Local beginning step position "
169 << startPointinLocal.x() << " " << startPointinLocal.y() << " "
170 << startPointinLocal.z());
171 ATH_MSG_DEBUG(" Local end step position "
172 << endPointinLocal.x() << " " << endPointinLocal.y() << " "
173 << endPointinLocal.z());
174#endif
175
176 // compute number of sub steps
177 // = 1 if no charge collection
178 // otherwise 200 microns (dstep) division
179
180 G4int nsub_step=1;
181 if (m_IflCur) nsub_step=(int) (thisStepLength*(1./dstep)) + 1;
182 // delta is fraction of step between two sub steps
183 G4double delta=1./((double) nsub_step);
184#ifdef DEBUGSTEP
185 ATH_MSG_DEBUG(" nsub_step,delta " << nsub_step << " " << delta);
186#endif
187
188 G4double energy = a_step->GetTotalEnergyDeposit() * a_step->GetTrack()->GetWeight();
189 if (m_BirksLaw) {
190 const double EField = 10. ; // 10 kV/cm electric field in presampler gap
191 const double wholeStepLengthCm = a_step->GetStepLength() / CLHEP::cm;
192 energy = (*m_birksLaw)(energy, wholeStepLengthCm, EField);
193 }
194
195 // loop over sub steps
196
197
198 // share total step energy evenly in each sub step
199 energy = energy / ((float) (nsub_step));
200 for (G4int i=0;i<nsub_step;i++) {
201
202#ifdef DEBUGSTEP
203 ATH_MSG_DEBUG(" Energy for sub step " << energy);
204#endif
205
206 // position for this substep
207 G4double fraction=(((G4double) i)+0.5)*delta;
208 G4ThreeVector subinLocal=startPointinLocal*(1-fraction) + endPointinLocal*fraction;
209 G4double xloc = subinLocal.x();
210 G4double yloc = subinLocal.y();
211 G4double zloc = subinLocal.z();
212 // call geometry method to find cell from local position
213 // status = true if hit is in normal region (13mm LAr gap)
214 // this method fills the cell number as well as coordinates in the electrode frame
215 LArG4::BarrelPresampler::CalcData currentCellData;
216 bool status = m_geometry->findCell(currentCellData,xloc,yloc,zloc);
217
218 // check fiducical cuts
219 if (status) {
220
221 // compute cell identifier
222 G4int zSide;
223 if (m_testbeam)
224 zSide = 1;
225 else
226 zSide = ( startPoint.z() > 0.) ? 1 : -1;
227
228 G4int region = currentCellData.region;
229 G4int etaBin = currentCellData.etaBin;
230 G4int phiBin = currentCellData.phiBin;
231 G4int sampling = currentCellData.sampling;
232
233 if( zSide == -1 )
234 {
235 // following code for an Y-axis rotation to define Z < 0. Barrel part
236 phiBin = 31 - phiBin;
237 if(phiBin < 0 ) phiBin += 64;
238 }
239
240 // check identifier
241 if (sampling !=0 || region != 0 ||
242 etaBin <0 || etaBin > 60 || phiBin <0 || phiBin>63) continue;
243
244 // fill identifier
245 identifier2.clear();
246 identifier2 << 4 // LArCalorimeter
247 << 1 // LArEM
248 << zSide
249 << sampling
250 << region
251 << etaBin
252 << phiBin;
253
254 // time computation is not necessarily correct for test beam
255 G4double time;
256 if (m_testbeam)
257 {
258 time=0.;
259 }
260 else
261 {
262 G4double tof;
263 tof = thisStepPoint->GetGlobalTime() / Units::ns;
264 time = tof - thisStepPoint->GetPosition().mag() * (1. / (CLHEP::c_light * CLHEP::ns));
265 }
266
267 G4double Current;
268 if (!m_IflCur) {
269 // no charge collection Current=E from Geant
270 Current=energy;
271 }
272 else {
273 // get module number and coordinates in electrode frame
274 G4int imodule = currentCellData.module;
275 G4double x0 = currentCellData.distElec;
276 G4double y0 = currentCellData.xElec;
277 // full symmetry for angle=0 electrodes
278 if (imodule>1) {
279 x0=std::fabs(x0);
280 y0=std::fabs(y0);
281 }
282 // reduced symmetry (point symmetry around 0) for tilted electrodes
283 if (imodule==0 || imodule ==1) {
284 if (x0<0) {
285 x0=-1.*x0;
286 y0=-1.*y0;
287 }
288 }
289#ifdef DEBUGSTEP
290 ATH_MSG_DEBUG(" set current map for module " << imodule);
291#endif
292 const CurrMap* cm = m_psmap->GetMap(imodule);
293 if (!cm) {
294 ATH_MSG_WARNING(" LArBarrelPresamplerCalculator: cannot get map for module " << imodule);
295 continue;
296 }
297 double current0,current1,current2,gap;
298
299 // get information from current map
300 cm->GetAll(x0,y0,&gap,&current0,&current1,&current2);
301#ifdef DEBUGSTEP
302 ATH_MSG_DEBUG(" module,x0,y0,current0 from map " << imodule << " " << x0 << " " << y0 << " " << current0);
303#endif
304
305 // assume HV=2000 everywhere for the time being
306 Current = energy*current0;
307
308 }
309
310 // check if we have a new hit in a different cell, or if we add this substep
311 // to an already existing hit
312 G4bool found=false;
313 for (unsigned int j=0; j<hdata.size(); j++) {
314 if (hdata[j].id==identifier2) {
315 hdata[j].energy += Current;
316 hdata[j].time += time*Current;
317 found=true;
318 break;
319 }
320 } // loop over hits
321 if (!found) {
322 hdata.emplace_back(LArHitData{identifier2, time*Current, Current});
323 } // hit was not existing before
324
325
326 } // hit fulfills fiducial cuts
327 } // end loop over sub steps
328
329#ifdef DEBUGSTEP
330 ATH_MSG_DEBUG("Number of hits for this step " << m_nhits << " " << hdata.size());
331#endif
332
333 // finalise time computations
334 for (unsigned int i=0;i<hdata.size();i++) {
335 if (std::fabs(hdata[i].energy)>1e-6) hdata[i].time=hdata[i].time/hdata[i].energy;
336 else hdata[i].time=0.;
337#ifdef DEBUGSTEP
338 ATH_MSG_DEBUG("Hit Energy/Time " << hdata[i].energy << " " << hdata[i].time);
339#endif
340 }
341
342 if (!hdata.empty()) return true;
343 else return false;
344
345}
#define ATH_MSG_WARNING(x)
gap(flags, cells_name, *args, **kw)
time(flags, cells_name, *args, **kw)
status
Definition merge.py:16
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setPhiMap phiBin
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setInterceptInner setEtaMap etaBin

Member Data Documentation

◆ m_Birksk

Gaudi::Property<double> LArCalculatorSvcImp::m_Birksk {this, "Birksk", 0.05832}
protectedinherited

Definition at line 27 of file LArCalculatorSvcImp.h.

27{this, "Birksk", 0.05832};

◆ m_BirksLaw

Gaudi::Property<bool> LArCalculatorSvcImp::m_BirksLaw {this, "BirksLaw", true}
protectedinherited

Definition at line 23 of file LArCalculatorSvcImp.h.

23{this, "BirksLaw", true};

◆ m_birksLaw

const LArG4BirksLaw* LArBarrelPresamplerCalculator::m_birksLaw {nullptr}
private

Definition at line 61 of file LArBarrelPresamplerCalculator.h.

61{nullptr};

◆ m_detectorName

Gaudi::Property<std::string> LArBarrelPresamplerCalculator::m_detectorName {this, "DetectorName", "LArMgr"}
private

Definition at line 56 of file LArBarrelPresamplerCalculator.h.

56{this, "DetectorName", "LArMgr"};

◆ m_geometry

ServiceHandle<ILArBarrelPresamplerGeometry> LArBarrelPresamplerCalculator::m_geometry {this, "GeometryCalculator", "LArBarrelPresamplerGeometry"}
private

Definition at line 51 of file LArBarrelPresamplerCalculator.h.

51{this, "GeometryCalculator", "LArBarrelPresamplerGeometry"};

◆ m_IflCur

Gaudi::Property<bool> LArBarrelPresamplerCalculator::m_IflCur {this, "IflCur", true}
private

Definition at line 53 of file LArBarrelPresamplerCalculator.h.

53{this, "IflCur", true};

◆ m_OOTcut

Gaudi::Property<double> LArCalculatorSvcImp::m_OOTcut {this, "OOTcut", 300*CLHEP::ns}
protectedinherited

Definition at line 30 of file LArCalculatorSvcImp.h.

30{this, "OOTcut", 300*CLHEP::ns};

◆ m_psmap

const PsMap* LArBarrelPresamplerCalculator::m_psmap {nullptr}
private

Definition at line 60 of file LArBarrelPresamplerCalculator.h.

60{nullptr};

◆ m_testbeam

Gaudi::Property<bool> LArBarrelPresamplerCalculator::m_testbeam {this, "isTestbeam", false}
private

Definition at line 58 of file LArBarrelPresamplerCalculator.h.

58{this, "isTestbeam", false};

◆ m_volname

G4String LArBarrelPresamplerCalculator::m_volname {"LArMgr::LAr::Barrel::Presampler"}
private

Definition at line 62 of file LArBarrelPresamplerCalculator.h.

62{"LArMgr::LAr::Barrel::Presampler"};

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