ATLAS Offline Software
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes | List of all members
LArBarrelPresamplerCalculator Class Reference

#include <LArBarrelPresamplerCalculator.h>

Inheritance diagram for LArBarrelPresamplerCalculator:
Collaboration diagram for LArBarrelPresamplerCalculator:

Public Member Functions

 LArBarrelPresamplerCalculator (const std::string &name, ISvcLocator *pSvcLocator)
 
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

bool m_BirksLaw
 
double m_Birksk
 
double m_OOTcut
 

Private Member Functions

 LArBarrelPresamplerCalculator (const LArBarrelPresamplerCalculator &)
 
LArBarrelPresamplerCalculatoroperator= (const LArBarrelPresamplerCalculator &)
 

Private Attributes

ServiceHandle< ILArBarrelPresamplerGeometrym_geometry
 
const PsMapm_psmap
 
bool m_IflCur
 
const LArG4BirksLawm_birksLaw
 
std::string m_detectorName
 
bool m_testbeam
 
G4String m_volname
 

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  , m_geometry("LArBarrelPresamplerGeometry", name)
44  , m_psmap(nullptr)
45  , m_IflCur(true)
46  , m_birksLaw(nullptr)
47  , m_detectorName("LArMgr")
48  , m_testbeam(false)
49  , m_volname("LArMgr::LAr::Barrel::Presampler")
50 {
51  declareProperty("GeometryCalculator", m_geometry);
52  declareProperty("IflCur",m_IflCur);
53  declareProperty("DetectorName",m_detectorName);
54  declareProperty("isTestbeam",m_testbeam);
55 }

◆ LArBarrelPresamplerCalculator() [2/2]

LArBarrelPresamplerCalculator::LArBarrelPresamplerCalculator ( const LArBarrelPresamplerCalculator )
private

Member Function Documentation

◆ finalize()

StatusCode LArBarrelPresamplerCalculator::finalize ( )
finaloverridevirtual

Definition at line 98 of file LArBarrelPresamplerCalculator.cxx.

99 {
100  if (m_BirksLaw) delete m_birksLaw;
101  return StatusCode::SUCCESS;
102 }

◆ initialize()

StatusCode LArBarrelPresamplerCalculator::initialize ( )
finaloverridevirtual

Definition at line 57 of file LArBarrelPresamplerCalculator.cxx.

58 {
59  ATH_MSG_DEBUG("LArBarrelPresamplerCalculator: Beginning initialization ");
60  // Initialize private members.
61 
62  // Access source of detector parameters.
64 
65  // Initialize the geometry calculator
66  ATH_CHECK(m_geometry.retrieve());
67 
68  if (m_IflCur)
69  {
70  ATH_MSG_DEBUG(" LArBarrelPresamplerCalculator: start reading current maps");
72  }
73 
74  // Get the out-of-time cut from the detector parameters routine.
75  m_OOTcut = parameters->GetValue("LArExpHallOutOfTimeCut"); //FIXME should be done by configurables
76 
77  if (m_BirksLaw)
78  {
79  const double Birks_LAr_density = 1.396;
80  m_birksLaw = new LArG4BirksLaw(Birks_LAr_density,m_Birksk);
81  ATH_MSG_DEBUG(" LArBarrelPresamplerCalculator: Birks' law ON ");
82  ATH_MSG_DEBUG(" LArBarrelPresamplerCalculator: parameter k " << m_birksLaw->k());
83  }
84  else
85  {
86  ATH_MSG_DEBUG(" LArBarrelPresamplerCalculator: Birks' law OFF");
87  }
88 
89  if(m_detectorName.empty()) m_volname="LAr::Barrel::Presampler";
90  else m_volname=m_detectorName+"::LAr::Barrel::Presampler";
91 
92  ATH_MSG_DEBUG("End of LArBarrelPresamplerCalculator initialization ");
93 
94  return StatusCode::SUCCESS;
95 }

◆ initializeForSDCreation()

virtual void LArCalculatorSvcImp::initializeForSDCreation ( )
inlineoverridevirtualinherited

Reimplemented in LArBarrelCalculator, and LArHECWheelCalculator.

Definition at line 18 of file LArCalculatorSvcImp.h.

18 {};

◆ isInTime()

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

Definition at line 40 of file LArBarrelPresamplerCalculator.h.

41  {
42  return !(std::fabs(hitTime) > m_OOTcut);
43  }

◆ OOTcut()

virtual G4float LArBarrelPresamplerCalculator::OOTcut ( ) const
inlinefinaloverridevirtual

Definition at line 35 of file LArBarrelPresamplerCalculator.h.

35 { return m_OOTcut; }

◆ operator=()

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

◆ Process()

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

Definition at line 105 of file LArBarrelPresamplerCalculator.cxx.

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

Member Data Documentation

◆ m_Birksk

double LArCalculatorSvcImp::m_Birksk
protectedinherited

Definition at line 25 of file LArCalculatorSvcImp.h.

◆ m_BirksLaw

bool LArCalculatorSvcImp::m_BirksLaw
protectedinherited

Definition at line 22 of file LArCalculatorSvcImp.h.

◆ m_birksLaw

const LArG4BirksLaw* LArBarrelPresamplerCalculator::m_birksLaw
private

Definition at line 58 of file LArBarrelPresamplerCalculator.h.

◆ m_detectorName

std::string LArBarrelPresamplerCalculator::m_detectorName
private

Definition at line 61 of file LArBarrelPresamplerCalculator.h.

◆ m_geometry

ServiceHandle<ILArBarrelPresamplerGeometry> LArBarrelPresamplerCalculator::m_geometry
private

Definition at line 52 of file LArBarrelPresamplerCalculator.h.

◆ m_IflCur

bool LArBarrelPresamplerCalculator::m_IflCur
private

Definition at line 56 of file LArBarrelPresamplerCalculator.h.

◆ m_OOTcut

double LArCalculatorSvcImp::m_OOTcut
protectedinherited

Definition at line 28 of file LArCalculatorSvcImp.h.

◆ m_psmap

const PsMap* LArBarrelPresamplerCalculator::m_psmap
private

Definition at line 54 of file LArBarrelPresamplerCalculator.h.

◆ m_testbeam

bool LArBarrelPresamplerCalculator::m_testbeam
private

Definition at line 63 of file LArBarrelPresamplerCalculator.h.

◆ m_volname

G4String LArBarrelPresamplerCalculator::m_volname
private

Definition at line 65 of file LArBarrelPresamplerCalculator.h.


The documentation for this class was generated from the following files:
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
LArBarrelPresamplerCalculator::m_IflCur
bool m_IflCur
Definition: LArBarrelPresamplerCalculator.h:56
LArG4Identifier
Definition: LArG4Identifier.h:121
LArG4::BarrelPresampler::CalcData::xElec
G4double xElec
Definition: ILArBarrelPresamplerGeometry.h:34
LArG4::BarrelPresampler::CalcData::region
G4int region
Definition: ILArBarrelPresamplerGeometry.h:28
LArBarrelPresamplerCalculator::m_volname
G4String m_volname
Definition: LArBarrelPresamplerCalculator.h:65
LArG4::BarrelPresampler::CalcData
Definition: ILArBarrelPresamplerGeometry.h:26
LArGeo::VDetectorParameters
Definition: VDetectorParameters.h:29
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
LArG4::BarrelPresampler::CalcData::phiBin
G4int phiBin
Definition: ILArBarrelPresamplerGeometry.h:30
LArG4Identifier::clear
void clear()
LArHitData
Definition: ILArCalculatorSvc.h:23
LArCalculatorSvcImp::LArCalculatorSvcImp
LArCalculatorSvcImp(const std::string &name, ISvcLocator *pSvcLocator)
Definition: LArCalculatorSvcImp.cxx:9
LArG4::BarrelPresampler::CalcData::module
G4int module
Definition: ILArBarrelPresamplerGeometry.h:32
LArCalculatorSvcImp::m_BirksLaw
bool m_BirksLaw
Definition: LArCalculatorSvcImp.h:18
LArCalculatorSvcImp::m_Birksk
double m_Birksk
Definition: LArCalculatorSvcImp.h:25
CaloSwCorrections.gap
def gap(flags, cells_name, *args, **kw)
Definition: CaloSwCorrections.py:212
LArG4BirksLaw::k
double k() const
Definition: LArG4BirksLaw.h:19
cm
const double cm
Definition: Simulation/ISF/ISF_FastCaloSim/ISF_FastCaloSimParametrization/tools/FCAL_ChannelMap.cxx:25
LArG4::BarrelPresampler::CalcData::distElec
G4double distElec
Definition: ILArBarrelPresamplerGeometry.h:33
PsMap::GetMap
const CurrMap * GetMap(int module) const
Definition: PsMap.cxx:55
xAOD::etaBin
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setInterceptInner setEtaMap etaBin
Definition: L2StandAloneMuon_v1.cxx:148
ParticleGun_FastCalo_ChargeFlip_Config.energy
energy
Definition: ParticleGun_FastCalo_ChargeFlip_Config.py:78
CurrMap
Definition: CurrMap.h:10
lumiFormat.i
int i
Definition: lumiFormat.py:85
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
LArBarrelPresamplerCalculator::m_psmap
const PsMap * m_psmap
Definition: LArBarrelPresamplerCalculator.h:54
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
xAOD::double
double
Definition: CompositeParticle_v1.cxx:159
MuonR4::SegmentFit::ParamDefs::x0
@ x0
MuonR4::SegmentFit::ParamDefs::y0
@ y0
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:221
xAOD::phiBin
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setPhiMap phiBin
Definition: L2StandAloneMuon_v2.cxx:144
LArBarrelPresamplerCalculator::m_detectorName
std::string m_detectorName
Definition: LArBarrelPresamplerCalculator.h:61
python.PhysicalConstants.c_light
float c_light
Definition: PhysicalConstants.py:63
LArBarrelPresamplerCalculator::m_geometry
ServiceHandle< ILArBarrelPresamplerGeometry > m_geometry
Definition: LArBarrelPresamplerCalculator.h:52
python.SystemOfUnits.mm
int mm
Definition: SystemOfUnits.py:83
LArG4BirksLaw
Definition: LArG4BirksLaw.h:8
LArG4::BarrelPresampler::CalcData::sampling
G4int sampling
Definition: ILArBarrelPresamplerGeometry.h:27
CondAlgsOpts.found
int found
Definition: CondAlgsOpts.py:101
CaloSwCorrections.time
def time(flags, cells_name, *args, **kw)
Definition: CaloSwCorrections.py:242
LArG4::BarrelPresampler::CalcData::etaBin
G4int etaBin
Definition: ILArBarrelPresamplerGeometry.h:29
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
LArCalculatorSvcImp::m_OOTcut
double m_OOTcut
Definition: LArCalculatorSvcImp.h:28
LArBarrelPresamplerCalculator::m_testbeam
bool m_testbeam
Definition: LArBarrelPresamplerCalculator.h:63
physics_parameters.parameters
parameters
Definition: physics_parameters.py:144
python.SystemOfUnits.ns
int ns
Definition: SystemOfUnits.py:130
PsMap::GetPsMap
static const PsMap * GetPsMap()
Definition: PsMap.cxx:40
merge.status
status
Definition: merge.py:17
hitTime
float hitTime(const AFP_SIDSimHit &hit)
Definition: AFP_SIDSimHit.h:39
LArGeo::VDetectorParameters::GetInstance
static const VDetectorParameters * GetInstance()
Definition: VDetectorParameters.cxx:29
readCCLHist.float
float
Definition: readCCLHist.py:83
LArBarrelPresamplerCalculator::m_birksLaw
const LArG4BirksLaw * m_birksLaw
Definition: LArBarrelPresamplerCalculator.h:58