ATLAS Offline Software
BarrelFastSimDedicatedSD.cxx
Go to the documentation of this file.
1 /*
2  Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3 */
4 
6 #include "LArSimEvent/LArHit.h"
12 #include "LArG4Code/EnergySpot.h"
13 #include "StoreGate/StoreGateSvc.h"
14 #include "CLHEP/Geometry/Point3D.h"
15 #include "CLHEP/Geometry/Transform3D.h"
17 
18 using HepGeom::Point3D;
20 using CLHEP::Hep3Vector;
21 
22 // Constructor:
24  : IFastSimDedicatedSD("BarrelFastSimDedicatedSD", detStore)
25  , m_embManager(nullptr)
26  , m_accordionDetails(nullptr)
27  , m_absorberSections(nullptr)
28 {
29  G4cout << GetName() << "::initialize()" << G4endl;
30  if ( detStore->retrieve( m_embManager ).isFailure() ){
31  throw std::runtime_error("Could not retrieve EMB manager");
32  }
35 }
36 
37 // ProcessHitsMethod
39 
40  // Fill the identifier.
41  Point3D<double> globalPosition=spot.GetPosition();
42 
45 
46  Point3D<double> localPosition = globalPosition.z()<0 ? xfNeg*globalPosition : xfPos*globalPosition;
47  int zIndex = globalPosition.z() <0 ? 0:1;
48  double eta = localPosition.getEta();
49  double phi = localPosition.getPhi();
50  double r = localPosition.perp();
51 
52  bool implementWaves=true;
53 
56  const EMBDetectorRegion *region = *e;
57  if (region->getEndcapIndex()!=zIndex) continue;
58  const EMBDetDescr *regionDescriptor=region->getDescriptor();
59  const CellBinning & etaBinning=regionDescriptor->getEtaBinning();
60  if (eta>etaBinning.getStart() && eta<etaBinning.getEnd()) {
61  unsigned int etaIndex = int((eta - etaBinning.getStart())/etaBinning.getDelta()) + etaBinning.getFirstDivisionNumber();
62  unsigned int phiIndex = 0;
63  EMBCellConstLink cellPtr = region->getEMBCell(etaIndex,phiIndex);
64  double rmax = cellPtr->getRLocal(EMBCell::BACK);
65  double rmin = cellPtr->getRLocal(EMBCell::FRONT);
66  if (r>rmin && r<rmax) {
67  const CellBinning & phiBinning=regionDescriptor->getPhiBinning();
68  double minPhi=std::min(phiBinning.getStart(), phiBinning.getEnd());
69  double maxPhi=std::max(phiBinning.getStart(), phiBinning.getEnd());
70 
71  while (phi<minPhi) phi += 2*M_PI;
72  while (phi>maxPhi) phi -= 2*M_PI;
73  unsigned int regionIndex = region->getRegionIndex();
74  unsigned int samplingIndex = region->getSamplingIndex();
75  unsigned int phiIndex = int((phi - phiBinning.getStart())/phiBinning.getDelta()) + phiBinning.getFirstDivisionNumber();
76 
77  // Now get the nominal cell and see if the indices need to be corrected for the accordion waves:
78  if (samplingIndex!=0) {
79  EMBCellConstLink cellPtr=region->getEMBCell(etaIndex,phiIndex);
80 
81  double delta=-2.0*M_PI/1024/2;
82 
83  double phiLocalUpper=cellPtr->getPhiLocalUpper()-delta;
84  double phiLocalLower=cellPtr->getPhiLocalLower()-delta;
85 
86  while (phiLocalLower<0) phiLocalLower += 2.0*M_PI;
87  while (phiLocalUpper<0) phiLocalUpper += 2.0*M_PI;
88  while (phiLocalLower>2*M_PI) phiLocalLower -= 2.0*M_PI;
89  while (phiLocalUpper>2*M_PI) phiLocalUpper -= 2.0*M_PI;
90 
91  int accordionIndexUpper=int(phiLocalUpper/(2*M_PI)*1024+0.25); if (accordionIndexUpper==1024) accordionIndexUpper=0;
92  int accordionIndexLower=int(phiLocalLower/(2*M_PI)*1024+0.25); if (accordionIndexLower==1024) accordionIndexLower=0;
93 
94  Point3D<double> A0 =
95  ( Point3D<double>(m_absorberSections->XCent(accordionIndexLower,1),m_absorberSections->YCent(accordionIndexLower,1),0)+
96  Point3D<double>(m_absorberSections->XCent(accordionIndexUpper,1),m_absorberSections->YCent(accordionIndexUpper,1),0))/2.0;
97 
98  Point3D<double> A1=
99  (Point3D<double>(m_absorberSections->XCent(accordionIndexLower,12),m_absorberSections->YCent(accordionIndexLower,12),0)+
100  Point3D<double>(m_absorberSections->XCent(accordionIndexUpper,12),m_absorberSections->YCent(accordionIndexUpper,12),0))/2.0-A0;
101 
102  Point3D<double> P0=Point3D<double>(localPosition.x(), localPosition.y(),0)-A0;
103  int stackIndex=int(A1.dot(P0)/A1.mag2()*22.0 + 3.0)/2 ;
104 
105  if (stackIndex<0 || stackIndex>13) {
106  G4cout << "Warning, bad stack index " << stackIndex << ' ' << rmin << ' ' << r << ' ' << A0.perp() << ' ' << A0 << ' ' << A1 << ' ' << P0 << G4endl;
107  if (implementWaves) return;
108  } else {
109  double xcent = m_absorberSections->XCent (accordionIndexLower,stackIndex);
110  double ycent = m_absorberSections->YCent (accordionIndexLower,stackIndex);
111  double cosU = m_absorberSections->Cosu (accordionIndexLower,stackIndex);
112  double sinU = m_absorberSections->Sinu (accordionIndexLower,stackIndex);
113  double halfLength = m_absorberSections->HalfLength(accordionIndexLower,stackIndex);
114 
115  Point3D<double> u(xcent-halfLength*cosU, ycent-halfLength*sinU,0);
116  Point3D<double> v(xcent+halfLength*cosU, ycent+halfLength*sinU,0);
117  Point3D<double> x=v-u,y=Hep3Vector(localPosition.x(),localPosition.y(),0)-u;
118  if ((zIndex==0 && x.cross(y).z()>0) || (zIndex==1&& x.cross(y).z()<0)) {
119  if (implementWaves) {
120  if (phiIndex==0) phiIndex=phiBinning.getNumDivisions()-1;
121  else phiIndex--;
122  }
123  } else {
124  double xcent = m_absorberSections->XCent (accordionIndexUpper,stackIndex);
125  double ycent = m_absorberSections->YCent (accordionIndexUpper,stackIndex);
126  double cosU = m_absorberSections->Cosu (accordionIndexUpper,stackIndex);
127  double sinU = m_absorberSections->Sinu (accordionIndexUpper,stackIndex);
128  double halfLength = m_absorberSections->HalfLength(accordionIndexUpper,stackIndex);
129 
130  Point3D<double> u(xcent-halfLength*cosU, ycent-halfLength*sinU,0);
131  Point3D<double> v(xcent+halfLength*cosU, ycent+halfLength*sinU,0);
132  Point3D<double> x=v-u,y=Hep3Vector(localPosition.x(),localPosition.y(),0)-u;
133  if ((zIndex==0 && x.cross(y).z()<0) || (zIndex==1 && x.cross(y).z()>0)) {
134  if (implementWaves) {
135  if (phiIndex==phiBinning.getNumDivisions()-1) phiIndex=0;
136  else phiIndex++;
137  }
138  }
139  }
140  }
141  }
142 
143  //static LArG4Identifier id;
144  m_larID.clear();
145  m_larID << 4 // LArCalorimeter
146  << 1 // LArEM
147  << ((zIndex==0) ? -1:1)
148  << samplingIndex
149  << regionIndex
150  << etaIndex
151  << phiIndex;
152  // call process to add this to the collection
153  SimpleHit(m_larID, spot.GetTime(), spot.GetEnergy());
154 
155  return;
156  }
157  }
158  }
159 
160 }
EMBDetectorRegion::getDescriptor
const EMBDetDescr * getDescriptor() const
Returns the Descriptor for this region.
Definition: EMBDetectorRegion.h:112
EMBDetDescr
Descriptor for regions of the electromagnetic barrel calorimeter.
Definition: EMBDetDescr.h:27
BarrelFastSimDedicatedSD::m_accordionDetails
const EMBAccordionDetails * m_accordionDetails
Definition: BarrelFastSimDedicatedSD.h:42
beamspotman.r
def r
Definition: beamspotman.py:676
max
#define max(a, b)
Definition: cfImp.cxx:41
phi
Scalar phi() const
phi method
Definition: AmgMatrixBasePlugin.h:64
EMBDetectorRegion.h
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
eta
Scalar eta() const
pseudorapidity method
Definition: AmgMatrixBasePlugin.h:79
EnergySpot
Definition: EnergySpot.h:18
GeoStraightAccSection::Sinu
const double & Sinu(int stackid, int cellid) const
BarrelFastSimDedicatedSD::m_embManager
const EMBDetectorManager * m_embManager
Definition: BarrelFastSimDedicatedSD.h:41
EMBDetectorRegion::getSamplingIndex
unsigned int getSamplingIndex() const
Returns the Sampling Layer Index.
Definition: EMBDetectorRegion.h:119
LArG4Identifier::clear
void clear()
M_PI
#define M_PI
Definition: ActiveFraction.h:11
IFastSimDedicatedSD
This is the interface for the fast simulation dedicated sensitive detector.
Definition: IFastSimDedicatedSD.h:13
BarrelFastSimDedicatedSD::m_absorberSections
const GeoStraightAccSection * m_absorberSections
Definition: BarrelFastSimDedicatedSD.h:43
x
#define x
EMBDetDescr::getPhiBinning
const CellBinning & getPhiBinning() const
The Binning in Phi.
Definition: EMBDetDescr.h:177
Trk::u
@ u
Enums for curvilinear frames.
Definition: ParamDefs.h:83
GeoStraightAccSection.h
GeoStraightAccSection::YCent
const double & YCent(int stackid, int cellid) const
StoreGateSvc
The Athena Transient Store API.
Definition: StoreGateSvc.h:128
BarrelFastSimDedicatedSD::BarrelFastSimDedicatedSD
BarrelFastSimDedicatedSD(StoreGateSvc *)
Definition: BarrelFastSimDedicatedSD.cxx:23
EnergySpot::GetEnergy
G4double GetEnergy() const
Definition: EnergySpot.h:37
BarrelFastSimDedicatedSD.h
Amg::Transform3D
Eigen::Affine3d Transform3D
Definition: GeoPrimitives.h:46
EMBDetectorRegion::getAbsoluteTransform
const Amg::Transform3D & getAbsoluteTransform(const GeoAlignmentStore *alignStore=nullptr) const
Returns the absolute transform of this element.
Definition: EMBDetectorRegion.cxx:27
EMBDetectorManager::getAccordionDetails
const EMBAccordionDetails * getAccordionDetails() const
Get accordion details class.
Definition: EMBDetectorManager.h:140
EnergySpot::GetTime
G4double GetTime() const
Definition: EnergySpot.h:43
EMBDetectorManager::beginDetectorRegion
EMBDetectorManager::DetectorRegionConstIterator beginDetectorRegion() const
Iterate over detector regions.
Definition: EMBDetectorManager.cxx:55
GeoStraightAccSection::Cosu
const double & Cosu(int stackid, int cellid) const
min
#define min(a, b)
Definition: cfImp.cxx:40
CLHEPtoEigenConverter.h
LArG4SimpleSD::SimpleHit
G4bool SimpleHit(const LArG4Identifier &lar_id, G4double time, G4double energy)
First method translates to this - also for fast sims.
fillTRTHists.zIndex
int zIndex
Definition: fillTRTHists.py:15
python.PyKernel.detStore
detStore
Definition: PyKernel.py:41
EMBCell::BACK
@ BACK
Definition: EMBCell.h:33
EMBDetDescr::getEtaBinning
const CellBinning & getEtaBinning() const
The Binning in Eta.
Definition: EMBDetDescr.h:184
EMBDetectorRegion::getEndcapIndex
EMBDetectorRegion::DetectorSide getEndcapIndex() const
The endcap index.
Definition: EMBDetectorRegion.h:163
EMBDetectorRegion::getEMBCell
EMBCellConstLink getEMBCell(unsigned int ieta, unsigned int iphi) const
Access to Cells.
Definition: EMBDetectorRegion.cxx:22
EMBDetectorRegion::getRegionIndex
unsigned int getRegionIndex() const
Returns the Region Index.
Definition: EMBDetectorRegion.h:126
eflowRec::phiIndex
unsigned int phiIndex(float phi, float binsize)
calculate phi index for a given phi
Definition: EtaPhiLUT.cxx:23
EnergySpot.h
python.PyAthena.v
v
Definition: PyAthena.py:157
BarrelFastSimDedicatedSD::ProcessSpot
void ProcessSpot(const EnergySpot &spot) override final
ProcessHitsMethod.
Definition: BarrelFastSimDedicatedSD.cxx:38
LArHit.h
BarrelFastSimDedicatedSD::m_larID
LArG4Identifier m_larID
My LAr identifier.
Definition: BarrelFastSimDedicatedSD.h:46
DiTauMassTools::MaxHistStrategyV2::e
e
Definition: PhysicsAnalysis/TauID/DiTauMassTools/DiTauMassTools/HelperFunctions.h:26
y
#define y
LArHitContainer.h
EMBDetectorManager::DetectorRegionConstIterator
std::vector< const EMBDetectorRegion * >::const_iterator DetectorRegionConstIterator
Definition: EMBDetectorManager.h:35
GeoStraightAccSection::HalfLength
const double & HalfLength(int stackid, int cellid) const
Amg::EigenTransformToCLHEP
HepGeom::Transform3D EigenTransformToCLHEP(const Amg::Transform3D &eigenTransf)
Converts an Eigen-based Amg::Transform3D into a CLHEP-based HepGeom::Transform3D.
Definition: CLHEPtoEigenConverter.h:120
EMBAccordionDetails::getAbsorberSections
const GeoStraightAccSection * getAbsorberSections() const
Absorber position details:
Definition: EMBAccordionDetails.cxx:298
GeoStraightAccSection::XCent
const double & XCent(int stackid, int cellid) const
EnergySpot::GetPosition
G4ThreeVector GetPosition() const
Definition: EnergySpot.h:40
EMBDetectorManager::getDetectorRegion
const EMBDetectorRegion * getDetectorRegion(unsigned int endcap, unsigned int sampling, unsigned int region) const
Random Access to detector regions.
Definition: EMBDetectorManager.cxx:65
EMBAccordionDetails.h
EMBDetectorRegion
Definition: EMBDetectorRegion.h:28
StoreGateSvc.h
EMBDetectorManager::endDetectorRegion
EMBDetectorManager::DetectorRegionConstIterator endDetectorRegion() const
Iterate over detector regions.
Definition: EMBDetectorManager.cxx:60
EMBDetectorManager.h
EMBCell::FRONT
@ FRONT
Definition: EMBCell.h:33