ATLAS Offline Software
Loading...
Searching...
No Matches
G4PolyconeGeoIDSvc.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
5// class header include
7
8// C includes
9#include <math.h>
10
11// framework includes
12#include "GaudiKernel/Bootstrap.h"
13#include "GaudiKernel/ISvcLocator.h"
14
15// STL includes
16#include <limits>
17#include <iostream>
18
19// ISF includes
21
22// Geant4 includes
23#include "geomdefs.hh"
24#include "G4Polycone.hh"
25#include "G4UnionSolid.hh"
26
27
28// DetectorDescription
31
33ISF::G4PolyconeGeoIDSvc::G4PolyconeGeoIDSvc(const std::string& name,ISvcLocator* svc) :
34 base_class(name,svc),
35 m_volume(),
37{
38}
39
40
45
46
47// Athena algtool's Hooks
49{
50 ATH_MSG_VERBOSE("initialize()");
51
52 // retrieve envelope definition service
53 ATH_CHECK( m_envDefSvc.retrieve() );
54
55 // create internal volume representations for the given dimensions
57 {
58 if ( createVolume( AtlasDetDescr::AtlasRegion(geoID)).isFailure())
59 {
60 ATH_MSG_ERROR("Unable to create volume representation for geoID="<<geoID);
61 }
62 }
63
64 // setup the converter for: G4 enum to ISF::InsideType
68
69 return StatusCode::SUCCESS;
70}
71
72
74{
75 // initialized with default return value
77
78 // check for volume to be present
79 G4VSolid *curVol = m_volume[geoID];
80 if (curVol)
81 {
82 G4ThreeVector posG4( pos.x(), pos.y(), pos.z() );
83 EInside g4Where = curVol->Inside( posG4 );
84
85 // use this trick to quickly convert G4 type to ISF::InsideType
86 where = m_typeConverter[g4Where];
87 }
88
89 // hand back the position information
90 return where;
91}
92
93
95{
96 // loop over geoIDs
97 for (unsigned short geoID=AtlasDetDescr::fFirstAtlasRegion; geoID<AtlasDetDescr::fNumAtlasRegions; ++geoID)
98 {
99 // point inside given geoID -> return geoID
101 {
102 return AtlasDetDescr::AtlasRegion(geoID);
103 }
104 }
105 // point is not inside any of the registered geoIDs
107}
108
109
111 const Amg::Vector3D &dir) const
112{
113
114 // if position inside a volume -> return that volume as the next geoID
116 if (validAtlasRegion(geoID))
117 {
118 return geoID;
119 }
120 //not inside - make a step towards the p direction and check if inside again
121 const Amg::Vector3D &dirUnit = dir.unit();
122 const Amg::Vector3D &posStep = pos+dirUnit;
123 AtlasDetDescr::AtlasRegion closestGeoID = identifyGeoID(posStep);
124 if (validAtlasRegion(closestGeoID))
125 {
126 return closestGeoID;
127 }
129}
130
131
133{
134
135 // ensure a proper numeric value for geoID
136 assert(validAtlasRegion(geoID));
137
138 ATH_MSG_INFO( "Building envelope volume for '" << AtlasDetDescr::AtlasRegionHelper::getName(geoID)
139 << "' (GeoID="<< geoID << ").");
140
141 // resetting volume for current geoID
142 m_volume[geoID] = 0;
143
144 // volume definitions
145 const double phimin = 0.;
146 const double deltaphi = 360.*CLHEP::deg;
147 const std::string volumeName( AtlasDetDescr::AtlasRegionHelper::getName(geoID));
148
149 {
150 // retrieve the vector of (r,z) values from the EnvelopeDefSvc
151 const RZPairVector &rz = m_envDefSvc->getRZBoundary( geoID );
152
153 // entries in the database table
154 size_t curVolNumPoints = rz.size();
155 if ( !curVolNumPoints)
156 {
157 ATH_MSG_ERROR("No entries for " << AtlasDetDescr::AtlasRegionHelper::getName(geoID) << " envelope in envelope definition service. Unable to build envelope.");
158 // resources cleanup
159 return StatusCode::FAILURE;
160 }
161
162 // arrays of (r,z) coordinates of the current envelope volume
163 double *z = new double[curVolNumPoints];
164 double *r = new double[curVolNumPoints];
165
166 // convert vector of pairs into two arrays
167 RZPairVector::const_iterator rzIt = rz.begin();
168 for ( size_t i=0; i<curVolNumPoints; ++i)
169 {
170 r[i] = rzIt->first;
171 z[i] = rzIt->second;
172
173 ++rzIt;
174 }
175
176 // create a G4Polycone volume with the given dimensions
177 std::stringstream curName;
178 curName << volumeName;
179 // store it in the m_volume array
180 m_volume[geoID] = new G4Polycone( curName.str(), phimin, deltaphi, curVolNumPoints, r, z);
181 // set the volume name
182 m_volume[geoID]->SetName( volumeName);
183
184 // cleanup
185 delete[] r;
186 delete[] z;
187 }
188
189 return StatusCode::SUCCESS;
190}
#define ATH_CHECK
Evaluate an expression and check for errors.
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
#define ATH_MSG_VERBOSE(x)
bool validAtlasRegion(AtlasDetDescr::AtlasRegion region)
Check a given AtlasRegion for its validity.
Definition AtlasRegion.h:40
static Double_t rz
std::vector< RZPair > RZPairVector
Definition RZPair.h:18
#define z
static const char * getName(int region)
ISF::InsideType m_typeConverter[ISF::fNumInsideTypes]
a quick way to convert G4 enum EInside to ISF::InsideType
AtlasDetDescr::AtlasRegion identifyGeoID(const Amg::Vector3D &pos) const
A static filter that returns the AtlasRegion of the given position.
ISF::InsideType inside(const Amg::Vector3D &pos, AtlasDetDescr::AtlasRegion geoID) const
Checks if the given position (ISFParticle) is inside a given AtlasRegion.
ServiceHandle< IEnvelopeDefSvc > m_envDefSvc
service providing the envelope dimensions for the different sub-detectors
G4VSolid * m_volume[AtlasDetDescr::fNumAtlasRegions]
G4PolyconeGeoIDSvc(const std::string &name, ISvcLocator *svc)
Constructor with parameters.
AtlasDetDescr::AtlasRegion identifyNextGeoID(const Amg::Vector3D &pos, const Amg::Vector3D &dir) const
Find the AtlasRegion that the particle will enter with its next infinitesimal step along the given di...
StatusCode createVolume(AtlasDetDescr::AtlasRegion geoID)
Retrieve and fill in the dimensions for the different AtlasRegion.
int r
Definition globals.cxx:22
Eigen::Matrix< double, 3, 1 > Vector3D
AtlasRegion
A simple enum of ATLAS regions and sub-detectors.
Definition AtlasRegion.h:21
InsideType
Definition IGeoIDSvc.h:22
@ fInside
Definition IGeoIDSvc.h:25
@ fOutside
Definition IGeoIDSvc.h:23
@ fSurface
Definition IGeoIDSvc.h:24