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