ATLAS Offline Software
LArBarrelCalculator.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 // The Cell Identifier for the EM Barrel readout cells
6 
7 // Adapted from code written by Gaston Parrour
8 // Adaptation by Sylvain Negroni
9 // 17-11-2003: G.Unal cleanup
10 // 18-03-2003: G.Unal major revision to use new current maps
11 
12 // #define DEBUGSTEP
13 // #define DEBUGSTEP2
14 
15 #include "LArBarrelCalculator.h"
16 #include "AccMap.h"
17 #include "MapEta.h"
18 
21 
22 //#include "LArG4RunControl/LArG4BarrelOptions.h"
23 
24 #include "G4ThreeVector.hh"
25 #include "G4StepPoint.hh"
26 #include "G4Step.hh"
27 #include "G4ios.hh"
28 #include "G4AffineTransform.hh"
29 #include "G4NavigationHistory.hh"
30 #include "G4VTouchable.hh"
31 #include "G4TouchableHistory.hh"
32 
33 #include <cmath>
34 #include <iostream>
35 
36 #include "GaudiKernel/ISvcLocator.h"
37 #include "GaudiKernel/Bootstrap.h"
38 #include "StoreGate/StoreGateSvc.h"
39 #include "AthenaKernel/Units.h"
40 
41 #include "LArHV/LArHVManager.h"
42 #include "LArHV/EMBHVManager.h"
43 #include "LArHV/EMBHVModule.h"
44 #include "LArHV/EMBHVElectrode.h"
45 
46 namespace Units = Athena::Units;
47 
48 // ================================================================================
49 LArBarrelCalculator::LArBarrelCalculator(const std::string& name, ISvcLocator* pSvcLocator)
50  : LArCalculatorSvcImp(name, pSvcLocator)
51  , m_geometry("LArBarrelGeometry", name)
52  , m_accmap(nullptr)
53  , m_IflCur(true)
54  , m_IflMapTrans(true)
55  , m_IflXtalk(true)
56  , m_dstep(.2*CLHEP::mm)
57  , m_birksLaw(nullptr)
58  , m_doHV(false)
59  , m_etaMaxBarrel(0)
60  , m_zMinBarrel(0)
61  , m_zMaxBarrel(0)
62  , m_rMinAccordion(0)
63  , m_rMaxAccordion(0)
64  , m_ThickAbs(0)
65  , m_ThickEle(0)
66  , m_NCellTot(0)
67  , m_NCellMax(0)
68  , m_testbeam(false)
69 {
70  ATH_MSG_DEBUG("LArBarrelCalculator: Beginning construction ");
71 
72  declareProperty("GeometryCalculator", m_geometry);
73  // define RUN conditions for the Barrel
74  declareProperty("EMBCurr",m_IflCur);
75  declareProperty("EMBEtaTrans",m_IflMapTrans);
76  declareProperty("EMBXtalk",m_IflXtalk);
77  declareProperty("EMBdstep", m_dstep);
78  declareProperty("EMBHVEnable",m_doHV);
79 }
80 
82 {
83  ATH_MSG_DEBUG("LArBarrelCalculator: Beginning initialization ");
84  if (m_BirksLaw) {
85  const double Birks_LAr_density = 1.396;
86  m_birksLaw = new LArG4BirksLaw(Birks_LAr_density,m_Birksk);
87  }
88 
89  // Access source of detector parameters.
91 
92  // Get the out-of-time cut from the detector parameters routine.
93  m_OOTcut = parameters->GetValue("LArExpHallOutOfTimeCut"); //FIXME should be done via configurables
94  ATH_MSG_DEBUG("**** OutOfTime cut " << m_OOTcut);
95 
96  // Main Barrel parameters
97  // All the UNITS are implicitly the GEANT4 ONES e.g. mm, rad, etc ...
98  m_etaMaxBarrel = parameters->GetValue("LArEMBMaxEtaAcceptance");
99  m_zMinBarrel = parameters->GetValue("LArEMBfiducialMothZmin");
100  m_zMaxBarrel = parameters->GetValue("LArEMBfiducialMothZmax");
101  m_NCellMax = (int) (parameters->GetValue("LArEMBnoOFPhysPhiCell"));
102  m_rMinAccordion = parameters->GetValue("LArEMBRadiusInnerAccordion");
103  m_rMaxAccordion = parameters->GetValue("LArEMBFiducialRmax");
104 
105  // absorbers and electrodes thickness
106  m_ThickAbs = 0.5*( parameters->GetValue("LArEMBThickAbsGlue")
107  +parameters->GetValue("LArEMBThickAbsIron")
108  +parameters->GetValue("LArEMBThickAbsLead"));
109 
110  G4double check = 0.5*( parameters->GetValue("LArEMBThinAbsGlue")
111  +parameters->GetValue("LArEMBThinAbsIron")
112  +parameters->GetValue("LArEMBThinAbsLead"));
113 
114  if (fabs(check-m_ThickAbs)>0.001)
115  {
116  ATH_MSG_WARNING("Thin and Thick Abs have difference thickness!");
117  }
118 
119  m_ThickEle= 0.5*( parameters->GetValue("LArEMBThickElecCopper")
120  +parameters->GetValue("LArEMBThickElecKapton"));
121 
122  // === GU 11/06/2003 total number of cells in phi
123  // to distinguish 1 module (testbeam case) from full Atlas
124  m_NCellTot = (int) (parameters->GetValue("LArEMBnoOFPhysPhiCell"));
125  if (m_NCellTot != 1024)
126  {
127  m_NCellMax=1024;
128  m_testbeam=true;
129  }
130  // ===
131  // access geometry computation class
132  ATH_CHECK(m_geometry.retrieve());
133 
134  // access current maps if required
135  if (m_IflCur)
136  {
137  ATH_MSG_DEBUG(" LArBarrelCalculator: start reading of current maps");
139  m_etamap1 = std::make_unique<MapEta>(1);
140  m_etamap2 = std::make_unique<MapEta>(2);
141  if (m_IflMapTrans) m_etamap3 = std::make_unique<MapEta>(3);
142  ATH_MSG_DEBUG(" LArBarrelCalculator: end of reading current maps");
143  }
144  // Initialize HV values
145  this->InitHV();
146 
147  ATH_MSG_DEBUG(" LArBarrelCalculator: s_NCellMax " << m_NCellMax);
148  ATH_MSG_DEBUG(" LArBarrelCalculator: s_NCellTot " << m_NCellTot);
149  ATH_MSG_DEBUG(" LArBarrelCalculator: s_rMinAccordion " << m_rMinAccordion);
150  ATH_MSG_DEBUG(" LArBarrelCalculator: s_rMaxAccordion " << m_rMaxAccordion);
151  ATH_MSG_DEBUG(" LArBarrelCalculator: s_zMinBarrel " << m_zMinBarrel);
152  ATH_MSG_DEBUG(" LArBarrelCalculator: s_zMaxBarrel " << m_zMaxBarrel);
153  ATH_MSG_DEBUG(" LArBarrelCalculator: s_etaMaxBarrel " << m_etaMaxBarrel);
154  ATH_MSG_DEBUG(" LArBarrelCalculator: s_ThickAbs " << m_ThickAbs);
155  ATH_MSG_DEBUG(" LArBarrelCalculator: s_ThickEle " << m_ThickEle);
156  if(m_IflCur) ATH_MSG_DEBUG(" LArBarrelCalculator: Deposited Energy dE/dX Corrected ==> CURRENT Option ON");
157  else ATH_MSG_DEBUG(" LArBarrelCalculator: Crude Deposited Energy dE/dX NO CURRENT option");
158  if (m_IflCur && m_IflMapTrans) ATH_MSG_DEBUG(" LArBarrelCalculator: Compute effect of E field around eta=0.8 ");
159  else ATH_MSG_DEBUG(" LArBarrelCalculator: Ignore effect of E field around eta=0.8 ");
160  if(m_BirksLaw)
161  {
162  ATH_MSG_DEBUG(" LArBarrelCalculator: Birks' law ON ");
163  ATH_MSG_DEBUG(" LArBarrelCalculator: parameter k " << m_birksLaw->k());
164  }
165  else
166  {
167  ATH_MSG_DEBUG(" LArBarrelCalculator: Birks' law OFF");
168  }
169 
170  return StatusCode::SUCCESS;
171 }
172 
173 // ============================================================================
175 {
176  if (m_BirksLaw) delete m_birksLaw;
177  return StatusCode::SUCCESS;
178 }
179 
180 // =============================================================================
181 G4bool LArBarrelCalculator::Process(const G4Step* step, std::vector<LArHitData>& hdata) const
182 {
183 
184  hdata.clear();
185 
186  LArG4Identifier identifier2;
187  LArG4Identifier identifier_xt1;
188  LArG4Identifier identifier_xt2;
189 
190  // check the Step content is non trivial
191  G4double thisStepEnergyDeposit = step->GetTotalEnergyDeposit() * step->GetTrack()->GetWeight();
192  G4double thisStepLength = step->GetStepLength() / Units::mm;
193 
194 #ifdef DEBUGSTEP
195  ATH_MSG_DEBUG("****** LArBarrelCalculator: Step energy,length "
196  << thisStepEnergyDeposit << " " << thisStepLength);
197 #endif
198 
199  if(thisStepEnergyDeposit <= 0. || thisStepLength <= 0.)
200  {
201 #ifdef DEBUGSTEP2
202  ATH_MSG_DEBUG(" Invalid hit trivial content");
203 #endif
204  return false;
205  }
206 
207  // Get Step geometrical parameters (first and end)
208  G4StepPoint *thisStepPoint = step->GetPreStepPoint();
209  G4StepPoint *thisStepBackPoint = step->GetPostStepPoint();
210  G4ThreeVector startPoint = thisStepPoint->GetPosition();
211  G4ThreeVector endPoint = thisStepBackPoint->GetPosition();
212 
213 #ifdef DEBUGSTEP
214  ATH_MSG_DEBUG(" Beginning step position "
215  << startPoint.x() << " " << startPoint.y() << " " << startPoint.z());
216 #endif
217 
218  // find zside from volume name
219  G4int zSide = 1;
220  if (!m_testbeam) {
221  const G4NavigationHistory* g4navigation = thisStepPoint->GetTouchable()->GetHistory();
222  G4int ndep = g4navigation->GetDepth();
223  for (G4int ii=0;ii<=ndep;ii++) {
224  G4VPhysicalVolume* v1 = g4navigation->GetVolume(ii);
225  G4String vname = v1->GetName();
226  if ( vname.find("NegPhysical") != std::string::npos) zSide=-1;
227  }
228  }
229 
230 
231  // BACK directly into the LOCAL half_Barrel Z > 0. PART (mother "stac_phys1")
232 
233  const G4AffineTransform transformation =
234  thisStepPoint->GetTouchable()->GetHistory()->GetTopTransform();
235 
236  G4ThreeVector startPointinLocal =
237  transformation.TransformPoint(startPoint);
238  G4ThreeVector endPointinLocal =
239  transformation.TransformPoint (endPoint);
240 
241 #ifdef DEBUGSTEP
242  ATH_MSG_DEBUG(" Local beginning step position "
243  << startPointinLocal.x() << " " << startPointinLocal.y() << " "
244  << startPointinLocal.z());
245  ATH_MSG_DEBUG(" Local end step position "
246  << endPointinLocal.x() << " " << endPointinLocal.y() << " "
247  << endPointinLocal.z());
248 #endif
249 
250  G4double energy = step->GetTotalEnergyDeposit() * step->GetTrack()->GetWeight(); // Despite the name, this is only ionization.
251 
252  if (m_BirksLaw) {
253  const double EField = 10.; // kV/cm, assume constant for now
254  const double wholeStepLengthCm = step->GetStepLength() / CLHEP::cm;
255  energy = (*m_birksLaw)(energy, wholeStepLengthCm, EField);
256  }
257 
258 
259  // compute number of sub steps
260  // = 1 if no charge collection
261  // otherwise 200 microns (dstep) division
262 
263  G4int nsub_step=1;
264  if (m_IflCur) nsub_step=(int) (thisStepLength/m_dstep) + 1;
265  // delta is fraction of step between two sub steps
266  G4double delta=1./((double) nsub_step);
267 #ifdef DEBUGSTEP
268  ATH_MSG_DEBUG(" nsub_step,delta " << nsub_step << " " << delta);
269 #endif
270 
271  energy /= ((float) (nsub_step)); // Divide into substeps.
272 
273  // loop over sub steps
274 
275  for (G4int i=0;i<nsub_step;i++) {
276 
277  bool Xtalk = false;
278 
279  G4double fraction=(((G4double) i)+0.5)*delta;
280  G4ThreeVector subinLocal=startPointinLocal*(1-fraction) + endPointinLocal*fraction;
281  G4double xloc = subinLocal.x();
282  G4double yloc = subinLocal.y();
283  G4double zloc = subinLocal.z();
284  G4double etaloc = subinLocal.pseudoRapidity();
285  G4double philoc = subinLocal.phi();
286  if(philoc<0.) philoc = philoc + 2.*M_PI;
287  G4double radloc = sqrt( xloc*xloc + yloc*yloc );
288 #ifdef DEBUGSTEP
289  ATH_MSG_DEBUG(" local position sub_step "
290  << xloc << " " << yloc << " " << zloc);
291 #endif
292 
293  // apply fiducial cuts
294  if (!FiducialCuts(radloc,zloc,etaloc)) {
295 #ifdef DEBUGSTEP
296  ATH_MSG_DEBUG("LArBarrelCalculator: outside fiducial cuts");
297 #endif
298  continue;
299  }
300  LArG4::Barrel::CalcData currentCellData;
301  m_geometry->findCell(currentCellData,xloc,yloc,zloc,radloc,etaloc,philoc,m_IflCur);
302 
303  if (currentCellData.cellID == 0)
304  {
305 #ifdef DEBUGSTEP
306  ATH_MSG_DEBUG("LArBarrelCalculator: Invalid hit CELLID == 0 ");
307  ATH_MSG_DEBUG("x,y,z local " << xloc << " " << yloc << " " << zloc);
308  ATH_MSG_DEBUG("radius " << radloc <<" etaloc "<< etaloc << " philoc " << philoc);
309 #endif
310  continue;
311  }
312  G4int region = currentCellData.region;
313  G4int etaBin = currentCellData.etaBin;
314  G4int phiBin = currentCellData.phiBin;
315  G4int sampling = currentCellData.sampling;
316 
317  if( zSide == -1 )
318  {
319  // following code for an Y-axis rotation to define Z < 0. Barrel part
320  if( sampling == 1 && region ==0 )
321  {
322  phiBin = 31 - phiBin;
323  if(phiBin < 0 ) phiBin += 64;
324  }
325  if( sampling == 1 && region ==1 )
326  {
327  phiBin = 127 - phiBin;
328  if(phiBin < 0 ) phiBin += 256;
329  }
330  if( sampling >= 2 )
331  {
332  phiBin = 127 - phiBin;
333  if(phiBin < 0 ) phiBin += 256;
334  }
335  }
336 #ifdef DEBUGSTEP
337  ATH_MSG_DEBUG(" region,side,sampling,eta,phi " << region << " " << zSide << " "
338  << sampling << " " << etaBin << " " << phiBin);
339  ATH_MSG_DEBUG(" distance to electrode,abs " << currentCellData.distElec << " "
340  << currentCellData.distAbs);
341  ATH_MSG_DEBUG(" local coordinates " << currentCellData.x0 << " " << currentCellData.y0);
342 #endif
343 
344  if (std::fabs(currentCellData.distElec)< m_ThickEle ||
345  std::fabs(currentCellData.distAbs) < m_ThickAbs) {
346 #ifdef DEBUGSTEP
347  ATH_MSG_DEBUG(" hit in absorber or electrode radius:" << radloc);
348 #endif
349  continue;
350  }
351 
352  identifier2.clear();
353  identifier2 << 4 // LArCalorimeter
354  << 1 // LArEM
355  << zSide
356  << sampling
357  << region
358  << etaBin
359  << phiBin;
360 
361  // time computation is not necessarily correct for test beam
362  // G4double time;
363  // if (m_testbeam)
364  // {
365  // time=0.;
366  // }
367  // else
368  // {
369  // const G4double tof = thisStepPoint->GetGlobalTime() / Units::ns;
370  // time = tof - ( thisStepPoint->GetPosition().mag()/CLHEP::c_light ) / Units::ns;
371  // }
372 
373  const G4double time = (m_testbeam)? 0.0 :
374  ( thisStepPoint->GetGlobalTime() - ( thisStepPoint->GetPosition().mag()/CLHEP::c_light ) ) / Units::ns;
375 
376 #ifdef DEBUGSTEP
377  ATH_MSG_DEBUG(" Energy for sub step " << energy);
378 #endif
379 
380  G4double Current;
381  G4double Current_xt1,Current_xt2;
382  Current_xt1=0;
383  Current_xt2=0;
384  if (!m_IflCur) {
385  // no charge collection Current=E from Geant
386  Current=energy;
387  }
388  else {
389  // full charge collection
390  G4double xmap,ymap;
391  G4int nfold = currentCellData.nfold;
392  G4double x0=currentCellData.x0;
393  G4double y0=currentCellData.y0;
394  if (x0<1500 || x0>1960 || y0>30 || y0<-30) {
395  ATH_MSG_WARNING("weird x0,y0 " << x0 << " " << y0);
396  }
397 #ifdef DEBUGSTEP
398  G4double rr = sqrt(x0*x0+y0*y0);
399  ATH_MSG_DEBUG(" radius,rad0 " << radloc << " " << rr);
400 #endif
401  const CurrMap* mapOfCurrent(nullptr);
402  bool UseFold=false;
403  // Are we close to a fold ? (fold 0 has some pathology)
404  if ((x0 > m_accmap->GetXmin(nfold) || nfold==0) &&
405  x0 < m_accmap->GetXmax(nfold) &&
406  y0 > m_accmap->GetYmin(nfold) &&
407  y0 < m_accmap->GetYmax(nfold) &&
408  (nfold != 13 || x0 < 1957.5) &&
409  (nfold>0 || x0 < 1504.6) ) {
410  xmap=x0;
411  ymap=y0;
412 #ifdef DEBUGSTEP
413  ATH_MSG_DEBUG(" Map for fold xmap,ymap " << nfold << " " << xmap << " " << ymap);
414 #endif
415  UseFold=true;
416  G4int sampMap = currentCellData.sampMap;
417  G4int etaMap = currentCellData.etaMap;
418  mapOfCurrent = m_accmap->GetMap(nfold,region,sampMap,etaMap);
419  // catch problem to find map
420  if (!mapOfCurrent) {
421  ATH_MSG_WARNING(" Problem to access map fold = " << nfold);
422  ATH_MSG_WARNING(" region,sampling,eta,fold " << region << " " << sampMap << " "
423  << etaMap << " " << nfold);
424  return false;
425  }
426  }
427  else {
428  G4int n;
429  G4int nstraight = currentCellData.nstraight;
430  if (nstraight%2==0) n=22;
431  else n=21;
432  mapOfCurrent = m_accmap->GetMap(n,region,sampling,etaBin);
433  // catch problem to find map
434  if (!mapOfCurrent) {
435  ATH_MSG_WARNING(" Problem to access map straight = " << nstraight);
436  return false;
437  }
438  xmap = currentCellData.xl;
439  ymap = currentCellData.distElec;
440  // special case for first straight section, which is shorter
441  if (nstraight==0) xmap = 0.5*(xmap+1.);
442 #ifdef DEBUGSTEP
443  ATH_MSG_DEBUG(" Map for straight xl,delec " << xmap << " " << ymap);
444 #endif
445  } // fold or straight
446  double gap;
447  double current0,current1,current2;
448  // get current for elementary charge
449  mapOfCurrent->GetAll(xmap,ymap,&gap,&current0,&current1,&current2);
450  G4double gap2=std::fabs(currentCellData.distElec)+std::fabs(currentCellData.distAbs)
452 
453  // in which HV cell are we ?
454  int ipm,ielec,ieta,iside;
455  if (zSide==1) ipm=1; // A side
456  else ipm=0; // C side
457  ielec=currentCellData.phiGap;
458  if (zSide==-1) {
459  ielec = 511 - ielec;
460  if(ielec < 0 ) ielec += 1024;
461  }
462  ieta=((int) (etaloc*(1./0.2)));
463  if (ieta>6) ieta=6; //part 1.4 to 1.475 is same HV as 1.2 to 1.4
464  iside=0; // phi lower than electrode 0, 1 for phi higher than electrode
465  if ((currentCellData.distElec>0 && zSide==1)
466  || (currentCellData.distElec<0 && zSide==-1) ) iside=1;
467 
468  // HV extrapolation
469  double current;
470  double hv=m_hv[ipm][ielec][ieta][iside];
471  if (hv>1995.) current=current0;
472  else if (hv>5.) current=ScaleHV(hv,current0,current1,current2);
473  else current=0.;
474 
475  // extrapolation for non nominal gap (allows to include sagging effect)
476  // i ~ (gap/gap2)**1.3
477  // gap = nominal gap from current map
478  // gap2 = real gap from geometry
479  // dgap = gap/gap2 -1
480  // at first order i ~ (1+1.3*dgap)
481  double dgap=0;
482  if (gap>1e-3 && gap2 >1e-3) dgap=gap/gap2-1;
483  current = current * (1. + 1.3*dgap);
484 
485 #ifdef DEBUGSTEP
486  ATH_MSG_DEBUG(" elementary current " << current0);
487  ATH_MSG_DEBUG(" Gap from Map/calculator " << gap << " " << gap2);
488 #endif
489 
490  Current = energy*current;
491 
492  // check if pathology...
493  if (!UseFold && std::fabs(currentCellData.distElec)>2.1 && current0 < 0.1) {
494  ATH_MSG_WARNING(" xl,distEle " << currentCellData.xl << " "
495  << currentCellData.distElec
496  << " str number " << currentCellData.nstraight
497  << " sampling,eta " << sampling << " " << etaBin << " "
498  << " current/E " << current0);
499  }
500 
501  // compute effect around transition at eta=0.8 from electric field
502  // read from map weighting factor ( E**1.3 = E*E**0.3)
503  // in x +-9mm around transition (assume symmetry around 0)
504  // in y distance from electrode in gap
505  bool InTrans=false;
506  //double Current_test=Current;
507  if (m_IflMapTrans) {
508  float etaTrans=0.8;
509  if (fabs(etaloc-etaTrans) < 0.025) {
510  double x=std::fabs(zloc-radloc*sinh(etaTrans))/cosh(etaTrans);
511  double y=std::fabs(currentCellData.distElec);
512  if ( x < m_etamap3->Xmax() ) {
513  double resp;
514  m_etamap3->GetData0(x,y,&resp);
515  Current = Current*resp;
516  InTrans=true;
517  }
518  } // eta = 0.8 +- 0.025
519  } // if m_IflMapTrans
520 
521  // simulate cross-talk effects and transitions in eta between cells
522  // (only in region 0, samplings 1 and 2)
523 
524  if (region==0) {
525 
526  double resp,xt0,xt1,xt2,deta;
527  if (sampling==1) {
528  deta=etaloc-0.003125*((double)etaBin+0.5);
529  m_etamap1->GetData(std::fabs(deta),std::fabs(currentCellData.distElec),
530  &resp,&xt0,&xt1,&xt2);
531 #ifdef DEBUGSTEP
532  ATH_MSG_DEBUG("hit in strip etaloc,etaBin,deta,delec,resp,xt0,xt1,xt2 "
533  << etaloc << " " << etaBin << " " << deta*1000 << " " << currentCellData.distElec
534  << " " << resp << " " << xt0 << " " << xt1 << " " << xt2);
535 #endif
536  if (!InTrans) Current = Current*resp;
537  if (m_IflXtalk && etaBin > 1 && etaBin < 446) {
538  Xtalk=true;
539  if (deta>0) {
540  identifier_xt1.clear();
541  identifier_xt1 << 4 << 1 << zSide << sampling << region << (etaBin+1) << phiBin;
542  Current_xt1 = Current*xt1;
543  identifier_xt2.clear();
544  identifier_xt2 << 4 << 1 << zSide << sampling << region << (etaBin-1) << phiBin;
545  Current_xt2 = Current*xt2;
546  }
547  else {
548  identifier_xt1.clear();
549  identifier_xt1 << 4 << 1 << zSide << sampling << region << (etaBin-1) << phiBin;
550  Current_xt1 = Current*xt1;
551  identifier_xt2.clear();
552  identifier_xt2 << 4 << 1 << zSide << sampling << region << (etaBin+1) << phiBin;
553  Current_xt2 = Current*xt2;
554  }
555  Current = Current*xt0;
556  } // m_IflXtalk = true
557  }
558  else if (sampling==2 && !InTrans) {
559  deta=etaloc-0.025*((double)etaBin+0.5);
560  m_etamap2->GetData0(std::fabs(deta),std::fabs(currentCellData.distElec),
561  &resp);
562 #ifdef DEBUGSTEP
563  ATH_MSG_DEBUG("hit in middle etaloc,etaBin,deta,delec,resp,xt0,xt1,xt2 "
564  << etaloc << " " << etaBin << " " << deta*1000 << " " << currentCellData.distElec
565  << " " << resp);
566 #endif
567  Current = Current*resp;
568  } // sampling =1 or 2
569 
570  } // region=0
571 
572  } // switch for current simulation
573 
574 
575  // check if we have a new hit in a different cell, or if we add this substep
576  // to an already existing hit
577  G4bool found=false;
578  for (unsigned int i=0; i<hdata.size(); i++) {
579  if (hdata[i].id==identifier2) {
580  hdata[i].energy += Current;
581  hdata[i].time += time*Current;
582  found=true;
583  break;
584  }
585  } // loop over hits
586  if (!found) {
587  LArHitData newdata = {identifier2, time*Current, Current};
588  hdata.push_back(newdata);
589  } // hit was not existing before
590 
591  if (Xtalk) {
592  for (unsigned int i=0; i<hdata.size(); i++) {
593  if (hdata[i].id==identifier_xt1) {
594  hdata[i].energy += Current_xt1;
595  hdata[i].time += time*Current_xt1;
596  found=true;
597  break;
598  }
599  } // loop over hits
600  if (!found) {
601  LArHitData newdata = {identifier_xt1, time*Current_xt1, Current_xt1};
602  hdata.push_back(newdata);
603  }
604  found=false;
605  for (unsigned int i=0; i<hdata.size(); i++) {
606  if (hdata[i].id==identifier_xt2) {
607  hdata[i].energy += Current_xt2;
608  hdata[i].time += time*Current_xt2;
609  found=true;
610  break;
611  }
612  } // loop over hits
613  if (!found) {
614  LArHitData newdata = {identifier_xt2, time*Current_xt2, Current_xt2};
615  hdata.push_back(newdata);
616  }
617  } // Xtalk true
618 
619 
620  } // *** End of loop over sub steps
621 
622 #ifdef DEBUGSTEP
623  ATH_MSG_DEBUG("Number of hits for this step " << m_nhits << " "
624  << hdata.size());
625 #endif
626 
627  // finalise time computations
628  for (unsigned int i=0;i<hdata.size();i++) {
629  if (std::fabs(hdata[i].energy)>1e-6) hdata[i].time=hdata[i].time/hdata[i].energy;
630  else hdata[i].time=0.;
631 #ifdef DEBUGSTEP
632  ATH_MSG_DEBUG("Hit Energy/Time "
633  << hdata[i].energy << " " << hdata[i].time);
634 #endif
635  }
636 
637  if (hdata.size()>0) return true;
638  return false;
639 }
640 
641 // ===============================================================================
642 // Cuts to define active region of accordion calorimeter
643 // 1500.024 < r < 1960.00 mm
644 // |eta| < 1.475
645 // 4 < z < 3164 mm (in local half barrel coordinates)
646 
647 G4bool LArBarrelCalculator::FiducialCuts(G4double radloc,G4double zloc,G4double etaloc) const
648 {
649  if (radloc < m_rMinAccordion ||
650  radloc > m_rMaxAccordion ||
651  etaloc > m_etaMaxBarrel ||
652  zloc < m_zMinBarrel ||
653  zloc > m_zMaxBarrel) return false;
654  else return true;
655 }
656 
657 // ==========================================================================
658 // HV values
659 // set all to 2000V. This should be changed if one wants to simulate
660 // some HV imperfections
662 {
663  ATH_MSG_DEBUG(" **** in LArBarrelCalculator::InitHV() ");
664 
665 
666 
667 
668  float defaultHvVal=2000.;
669  ATH_MSG_DEBUG(" defaultHvVal " << defaultHvVal);
670  for (int ipm=0;ipm<2;ipm++) {
671  for (int ielec=0;ielec<1024;ielec++) {
672  for (int ieta=0;ieta<7;ieta++) {
673  for (int iside=0;iside<2;iside++) {
674  m_hv[ipm][ielec][ieta][iside] = defaultHvVal;
675  }
676  }
677  }
678  }
679 
680  if (m_doHV) {
681  // initialize services
682  SmartIF<StoreGateSvc> pDetStore{Gaudi::svcLocator()->service("DetectorStore")};
683  if (!pDetStore) {
684  ATH_MSG_WARNING("LArBarrelCalculator::InitHV() unable to get Detector Store! Use default HV values.");
685  return;
686  }
687 
688  // get EMBHV Manager
689  const LArHVManager *manager = nullptr;
690  if (pDetStore->retrieve(manager)==StatusCode::SUCCESS) {
691  const EMBHVManager& hvManager=manager->getEMBHVManager();
692  const EMBHVManager::EMBHVData hvdata = hvManager.getDataSim();
693  ATH_MSG_DEBUG(" got HV Manager ");
694  // loop over HV modules
695  for (unsigned int iSide=0;iSide<2;iSide++) {
696  for (unsigned int iPhi=0;iPhi<16;iPhi++) {
697  for (unsigned int iSector=0;iSector<2;iSector++) {
698  for (unsigned int iEta=0;iEta<7;iEta++) {
699  const EMBHVModule& hvMod = hvManager.getHVModule(iSide,iEta+1,iPhi,iSector);
700  for (unsigned int ielec=0;ielec<32;ielec++) {
701  const EMBHVElectrode& electrode = hvMod.getElectrode(ielec);
702  unsigned jElec = ielec+32*iSector+64*iPhi;
703  for (unsigned int iGap=0;iGap<2;iGap++) {
704  double hv = hvdata.voltage (electrode, iGap);
705  ATH_MSG_DEBUG(" iSide,jElec,iEta,iGap,hv " << iSide << " " << jElec << " " << iEta << " " << iGap << " " << hv);
706  if (hv>-999.) m_hv[iSide][jElec][iEta][iGap] = hv;
707  }
708  }
709  }
710  }
711  }
712  }
713  }
714  else {
715  ATH_MSG_WARNING(" Unable to find HV Manager ");
716  }
717  }
718 
719  return;
720 }
721 
722 //==========================================================================
723 // extrapolate response with HV
724 // hv = actual HV
725 // curr0 = response for 2000V, curr1 = for 1000V, curr2 for 400V
726 // assumes resp = a* HV*b between 1000-2000 and 400-1000V
727 // (fixes b=0.57 for HV smaller than 400V)
728 
729 double LArBarrelCalculator::ScaleHV(double hv, double curr0, double curr1, double curr2) const
730 {
731  double b;
732  double resp=0.;
733  if (hv>1000.) {
734  if (std::fabs(curr1)>1e-6) {
735  double x=curr0/curr1;
736  if (x>1e-6) {
737  b=log(x)*(1./(log(2000.)-log(1000.)));
738  resp = curr0*exp(b*(log(hv)-log(2000.)));
739  }
740  }
741  }
742  else if (hv>400.) {
743  if (std::fabs(curr2)>1e-6) {
744  double x=curr1/curr2;
745  if (x>1e-6) {
746  b=log(x)*(1./(log(1000.)-log(400.)));
747  resp = curr1*exp(b*(log(hv)-log(1000.)));
748  }
749  }
750  }
751  else {
752  if (std::fabs(curr2)>1e-6) {
753  b=0.57;
754  resp = curr2*exp(b*(log(hv)-log(400.)));
755  }
756  }
757  return resp;
758 }
LArBarrelCalculator::Process
virtual G4bool Process(const G4Step *a_step, std::vector< LArHitData > &hdata) const override final
Definition: LArBarrelCalculator.cxx:181
AllowedVariables::e
e
Definition: AsgElectronSelectorTool.cxx:37
LArBarrelCalculator::m_geometry
ServiceHandle< ILArBarrelGeometry > m_geometry
Definition: LArBarrelCalculator.h:60
fillPileUpNoiseLumi.current
current
Definition: fillPileUpNoiseLumi.py:52
LArBarrelCalculator::m_etamap2
std::unique_ptr< MapEta > m_etamap2
Definition: LArBarrelCalculator.h:63
EMBHVElectrode
Definition: EMBHVElectrode.h:15
CurrMap::GetAll
void GetAll(double x, double y, double *gap, double *curr0, double *curr1, double *curr2) const
Definition: CurrMap.cxx:93
LArHVManager.h
LArG4BirksLaw.h
LArG4Identifier
Definition: LArG4Identifier.h:121
LArBarrelCalculator::m_ThickEle
double m_ThickEle
Definition: LArBarrelCalculator.h:85
LArBarrelCalculator::m_NCellTot
int m_NCellTot
Definition: LArBarrelCalculator.h:87
LArGeo::VDetectorParameters
Definition: VDetectorParameters.h:29
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
LArBarrelCalculator::m_accmap
const AccMap * m_accmap
Definition: LArBarrelCalculator.h:61
LArBarrelCalculator::m_etamap3
std::unique_ptr< MapEta > m_etamap3
Definition: LArBarrelCalculator.h:64
LArBarrelCalculator::m_etamap1
std::unique_ptr< MapEta > m_etamap1
Definition: LArBarrelCalculator.h:62
LArBarrelCalculator::m_IflCur
bool m_IflCur
Definition: LArBarrelCalculator.h:67
LArBarrelCalculator::m_etaMaxBarrel
double m_etaMaxBarrel
Definition: LArBarrelCalculator.h:77
LArBarrelCalculator::m_zMinBarrel
double m_zMinBarrel
Definition: LArBarrelCalculator.h:78
LArBarrelCalculator::m_hv
double m_hv[2][1024][7][2]
Definition: LArBarrelCalculator.h:100
LArG4Identifier::clear
void clear()
M_PI
#define M_PI
Definition: ActiveFraction.h:11
LArHitData
Definition: ILArCalculatorSvc.h:23
LArVG4DetectorParameters.h
LArBarrelCalculator::m_ThickAbs
double m_ThickAbs
Definition: LArBarrelCalculator.h:84
LArG4::Barrel::CalcData
Definition: ILArBarrelGeometry.h:26
LArBarrelCalculator::m_zMaxBarrel
double m_zMaxBarrel
Definition: LArBarrelCalculator.h:79
LArG4::Barrel::CalcData::phiGap
G4int phiGap
Definition: ILArBarrelGeometry.h:33
LArCalculatorSvcImp::m_BirksLaw
bool m_BirksLaw
Definition: LArCalculatorSvcImp.h:18
drawFromPickle.exp
exp
Definition: drawFromPickle.py:36
LArG4::Barrel::CalcData::x0
G4double x0
Definition: ILArBarrelGeometry.h:39
x
#define x
LArBarrelCalculator::InitHV
void InitHV()
Definition: LArBarrelCalculator.cxx:661
LArCalculatorSvcImp::m_Birksk
double m_Birksk
Definition: LArCalculatorSvcImp.h:25
LArG4::Barrel::CalcData::xl
G4double xl
Definition: ILArBarrelGeometry.h:38
CaloSwCorrections.gap
def gap(flags, cells_name, *args, **kw)
Definition: CaloSwCorrections.py:212
EMBHVElectrode.h
LArBarrelCalculator::m_NCellMax
int m_NCellMax
Definition: LArBarrelCalculator.h:88
LArBarrelCalculator::initialize
virtual StatusCode initialize() override final
Definition: LArBarrelCalculator.cxx:81
LArBarrelCalculator::m_testbeam
bool m_testbeam
Definition: LArBarrelCalculator.h:91
LArBarrelCalculator::LArBarrelCalculator
LArBarrelCalculator(const std::string &name, ISvcLocator *pSvcLocator)
Definition: LArBarrelCalculator.cxx:49
LArG4BirksLaw::k
double k() const
Definition: LArG4BirksLaw.h:19
EMBHVManager::getHVModule
const EMBHVModule & getHVModule(unsigned int iSide, unsigned int iEta, unsigned int iPhi, unsigned int iSector) const
Definition: EMBHVManager.cxx:194
cm
const double cm
Definition: Simulation/ISF/ISF_FastCaloSim/ISF_FastCaloSimParametrization/tools/FCAL_ChannelMap.cxx:25
LArG4::Barrel::CalcData::cellID
int cellID
Definition: ILArBarrelGeometry.h:27
LArBarrelCalculator::ScaleHV
double ScaleHV(double, double, double, double) const
Definition: LArBarrelCalculator.cxx:729
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
EMBHVManager::EMBHVData::voltage
double voltage(const EMBHVElectrode &electrode, const int &iGap) const
Definition: EMBHVManager.cxx:130
LArBarrelCalculator::m_birksLaw
const LArG4BirksLaw * m_birksLaw
Definition: LArBarrelCalculator.h:73
beamspotman.n
n
Definition: beamspotman.py:731
LArG4::Barrel::CalcData::sampling
G4int sampling
Definition: ILArBarrelGeometry.h:28
EL::StatusCode
::StatusCode StatusCode
StatusCode definition for legacy code.
Definition: PhysicsAnalysis/D3PDTools/EventLoop/EventLoop/StatusCode.h:22
ATH_MSG_DEBUG
#define ATH_MSG_DEBUG(x)
Definition: AthMsgStreamMacros.h:29
xAOD::etaMap
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setInterceptInner etaMap
Definition: L2StandAloneMuon_v1.cxx:144
CLHEP
STD'S.
Definition: IAtRndmGenSvc.h:19
LArBarrelCalculator::m_IflMapTrans
bool m_IflMapTrans
Definition: LArBarrelCalculator.h:68
LArG4::Barrel::CalcData::phiBin
G4int phiBin
Definition: ILArBarrelGeometry.h:31
LArBarrelCalculator::m_rMinAccordion
double m_rMinAccordion
Definition: LArBarrelCalculator.h:81
LArBarrelCalculator::m_IflXtalk
bool m_IflXtalk
Definition: LArBarrelCalculator.h:69
ATH_CHECK
#define ATH_CHECK
Definition: AthCheckMacros.h:40
LArG4::Barrel::CalcData::nfold
G4int nfold
Definition: ILArBarrelGeometry.h:35
AccMap.h
xAOD::double
double
Definition: CompositeParticle_v1.cxx:159
LArG4::Barrel::CalcData::etaMap
G4int etaMap
Definition: ILArBarrelGeometry.h:42
LArCalculatorSvcImp
Definition: LArCalculatorSvcImp.h:11
EMBHVManager::getDataSim
EMBHVData getDataSim() const
Definition: EMBHVManager.cxx:324
LArG4::Barrel::CalcData::sampMap
G4int sampMap
Definition: ILArBarrelGeometry.h:41
MapEta.h
LArBarrelCalculator::m_rMaxAccordion
double m_rMaxAccordion
Definition: LArBarrelCalculator.h:82
LArG4::Barrel::CalcData::region
G4int region
Definition: ILArBarrelGeometry.h:29
Athena::Units
Definition: Units.h:45
LArBarrelCalculator.h
EMBHVManager::EMBHVData
Definition: EMBHVManager.h:39
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:228
xAOD::phiBin
setSAddress setEtaMS setDirPhiMS setDirZMS setBarrelRadius setEndcapAlpha setEndcapRadius setPhiMap phiBin
Definition: L2StandAloneMuon_v2.cxx:144
plotBeamSpotMon.b
b
Definition: plotBeamSpotMon.py:77
LArBarrelCalculator::m_dstep
double m_dstep
Definition: LArBarrelCalculator.h:71
LArBarrelCalculator::FiducialCuts
G4bool FiducialCuts(G4double, G4double, G4double) const
Definition: LArBarrelCalculator.cxx:647
EMBHVManager
This class provides direct access to information on the HV electrodes within the barrels....
Definition: EMBHVManager.h:36
python.PhysicalConstants.c_light
float c_light
Definition: PhysicalConstants.py:63
Units.h
Wrapper to avoid constant divisions when using units.
AccMap::GetMap
const CurrMap * GetMap(int ifold, int ielecregion) const
Definition: AccMap.cxx:92
LArNewCalib_Delay_OFC_Cali.check
check
Definition: LArNewCalib_Delay_OFC_Cali.py:267
Trk::iPhi
@ iPhi
Definition: ParamDefs.h:47
LArG4::Barrel::CalcData::etaBin
G4int etaBin
Definition: ILArBarrelGeometry.h:30
LArG4::Barrel::CalcData::distElec
G4double distElec
Definition: ILArBarrelGeometry.h:36
python.SystemOfUnits.mm
int mm
Definition: SystemOfUnits.py:83
LArG4::Barrel::CalcData::nstraight
G4int nstraight
Definition: ILArBarrelGeometry.h:34
LArG4::Barrel::CalcData::y0
G4double y0
Definition: ILArBarrelGeometry.h:40
LArG4BirksLaw
Definition: LArG4BirksLaw.h:8
AccMap::GetAccMap
static const AccMap * GetAccMap()
Definition: AccMap.cxx:72
LArHVManager
This class provides access to the High Voltage throughout the LAr. High voltage conditions can also b...
Definition: LArHVManager.h:24
y
#define y
CondAlgsOpts.found
int found
Definition: CondAlgsOpts.py:101
CaloSwCorrections.time
def time(flags, cells_name, *args, **kw)
Definition: CaloSwCorrections.py:242
AccMap::GetYmin
float GetYmin(int ifold) const
Definition: AccMap.h:31
ATH_MSG_WARNING
#define ATH_MSG_WARNING(x)
Definition: AthMsgStreamMacros.h:32
LArCalculatorSvcImp::m_OOTcut
double m_OOTcut
Definition: LArCalculatorSvcImp.h:28
LArBarrelCalculator::m_doHV
bool m_doHV
Definition: LArBarrelCalculator.h:74
python.CaloCondTools.log
log
Definition: CaloCondTools.py:20
LArCellBinning.step
step
Definition: LArCellBinning.py:158
physics_parameters.parameters
parameters
Definition: physics_parameters.py:144
EMBHVModule::getElectrode
const EMBHVElectrode & getElectrode(unsigned int iElectrode) const
Definition: EMBHVModule.cxx:69
python.SystemOfUnits.ns
int ns
Definition: SystemOfUnits.py:130
python.Logging.manager
manager
Definition: PhysicsAnalysis/D3PDTools/AnaAlgorithm/python/Logging.py:92
rr
const boost::regex rr(r_r)
EMBHVModule
Describes one HV Module within the EMB.
Definition: EMBHVModule.h:20
LArG4::Barrel::CalcData::distAbs
G4double distAbs
Definition: ILArBarrelGeometry.h:37
xAOD::iEta
setScale setgFexType iEta
Definition: gFexJetRoI_v1.cxx:77
LArGeo::VDetectorParameters::GetInstance
static const VDetectorParameters * GetInstance()
Definition: VDetectorParameters.cxx:29
StoreGateSvc.h
AccMap::GetXmin
float GetXmin(int ifold) const
Definition: AccMap.h:25
readCCLHist.float
float
Definition: readCCLHist.py:83
EMBHVModule.h
EMBHVManager.h
LArBarrelCalculator::finalize
virtual StatusCode finalize() override final
Definition: LArBarrelCalculator.cxx:174