ATLAS Offline Software
Public Member Functions | Static Public Member Functions | List of all members
iGeant4::ISFG4GeoHelper Class Reference

#include <ISFG4GeoHelper.h>

Collaboration diagram for iGeant4::ISFG4GeoHelper:

Public Member Functions

 ISFG4GeoHelper ()=delete
 

Static Public Member Functions

static AtlasDetDescr::AtlasRegion nextGeoId (const G4Step *aStep, int truthVolLevel, ISF::IGeoIDSvc *geoIDSvc)
 
static bool checkVolumeDepth (G4LogicalVolume *logicalVol, int volLevel, int depth=0)
 
static AtlasDetDescr::AtlasRegion getNextGeoIDFromSvc (const G4StepPoint &postStep, const ISF::IGeoIDSvc &geoIDSvc)
 get the next GeoID using only the geoIDSvc More...
 

Detailed Description

Definition at line 24 of file ISFG4GeoHelper.h.

Constructor & Destructor Documentation

◆ ISFG4GeoHelper()

iGeant4::ISFG4GeoHelper::ISFG4GeoHelper ( )
delete

Member Function Documentation

◆ checkVolumeDepth()

bool iGeant4::ISFG4GeoHelper::checkVolumeDepth ( G4LogicalVolume *  logicalVol,
int  volLevel,
int  depth = 0 
)
static

Definition at line 109 of file ISFG4GeoHelper.cxx.

109  {
110 
111  //FIXME - can replace all this code with similar methods to those in MCTruthBase/src/RecordingEnvelope.cxx
112  if (lv==nullptr) { return false; }
113 
114  bool Cavern = false;
115  // Check the volumes rather explicitly
116 
117  if ( lv->GetName() == "BeamPipe::BeamPipe" ||
118  lv->GetName() == "IDET::IDET" ||
119  lv->GetName() == "ITK::ITK" ||
120  lv->GetName() == "CALO::CALO" ||
121  lv->GetName() == "MUONQ02::MUONQ02" ||
122  lv->GetName() == "TTR_BARREL::TTR_BARREL" ) {
123  if(depth!=volLevel) {
124  G4ExceptionDescription description;
125  description << "Volume " << lv->GetName() << " at depth " << depth << " instead of depth " << volLevel;
126  G4Exception("iGeant4::ISFG4GeoHelper", "checkVolumeDepth", FatalException, description);
127  return false; //The G4Exception call above should abort the job, but Coverity does not seem to pick this up.
128  }
129  //ATH_MSG_DEBUG("Volume " << lv->GetName() << " is correctly registered at depth " << depth);
130  }
131  else if ( lv->GetName()=="BeamPipe::BeamPipeCentral" ) { // Things that are supposed to be one deeper
132  if (depth!=volLevel+1) {
133  G4ExceptionDescription description;
134  description << "Volume " << lv->GetName() << " at depth " << depth << " instead of depth " << volLevel+1;
135  G4Exception("iGeant4::ISFG4GeoHelper", "checkVolumeDepth", FatalException, description);
136  return false; //The G4Exception call above should abort the job, but Coverity does not seem to pick this up.
137  }
138  //ATH_MSG_DEBUG("Volume " << lv->GetName() << " is correctly registered at depth " << depth);
139  }
140  else if ( lv->GetName().find("CavernInfra")!=std::string::npos ) { // Things that are supposed to be one shallower
141  if (depth==volLevel-1) {
142  Cavern=true;
143  //ATH_MSG_DEBUG("Volume " << lv->GetName() << " is correctly registered at depth " << depth);
144  // Note: a number of volumes exist with "CavernInfra" in the name at the wrong depth, so we just need to
145  // check that there's at least one at the right depth
146  } // Check of volume level
147  } // Check of volume name
148 
149  // Going through the volume depth
150  for (unsigned int i=0; i<lv->GetNoDaughters(); ++i) {
151  Cavern = Cavern || checkVolumeDepth( lv->GetDaughter(i)->GetLogicalVolume() , volLevel , depth+1 );
152  }
153  if (depth==0 && !Cavern && volLevel>1) {
154  G4ExceptionDescription description;
155  description << "No CavernInfra volume registered at depth " << volLevel-1;
156  G4Exception("iGeant4::ISFG4GeoHelper", "checkVolumeDepth", FatalException, description);
157  return false; //The G4Exception call above should abort the job, but Coverity does not seem to pick this up.
158  }
159 
160  return Cavern;
161 }

◆ getNextGeoIDFromSvc()

AtlasDetDescr::AtlasRegion iGeant4::ISFG4GeoHelper::getNextGeoIDFromSvc ( const G4StepPoint &  postStep,
const ISF::IGeoIDSvc geoIDSvc 
)
static

get the next GeoID using only the geoIDSvc

Definition at line 88 of file ISFG4GeoHelper.cxx.

89  {
90  const G4ThreeVector &postPos = postStep.GetPosition();
91  AtlasDetDescr::AtlasRegion nextGeoID = geoIDSvc.identifyGeoID(postPos.x(),
92  postPos.y(),
93  postPos.z());
94 
95  // if we ever run into problems with the current approach, the following
96  // takes the particle's traveling direction into account for finding the
97  // *next* volume it enters
98  //const G4ThreeVector &postMom = postStep->GetMomentum();
99  //nextGeoID = m_geoIDSvcQuick->identifyNextGeoID( postPos.x(),
100  // postPos.y(),
101  // postPos.z(),
102  // postMom.x(),
103  // postMom.y(),
104  // postMom.z() );
105 
106  return nextGeoID;
107 }

◆ nextGeoId()

AtlasDetDescr::AtlasRegion iGeant4::ISFG4GeoHelper::nextGeoId ( const G4Step *  aStep,
int  truthVolLevel,
ISF::IGeoIDSvc geoIDSvc 
)
static

Definition at line 20 of file ISFG4GeoHelper.cxx.

21 {
22  struct VolumeHolder {
23  G4LogicalVolume *BP, *ID, *CALO, *MU, *TTR;
24  };
25 
26  // one time initialization of volumes
27  static const VolumeHolder vh = [&]() {
28  VolumeHolder h{};
29  G4LogicalVolumeStore * lvs = G4LogicalVolumeStore::GetInstance();
30  for (size_t i=0;i<lvs->size();++i) {
31 
32  if ( !(*lvs)[i] ) { continue; }
33  if ( (*lvs)[i]->GetName() == "BeamPipe::BeamPipe" ) { h.BP = (*lvs)[i]; }
34  else if ( (*lvs)[i]->GetName() == "IDET::IDET" || (*lvs)[i]->GetName() == "ITK::ITK" ) { h.ID = (*lvs)[i]; }
35  else if ( (*lvs)[i]->GetName() == "CALO::CALO" ) { h.CALO = (*lvs)[i]; }
36  else if ( (*lvs)[i]->GetName() == "MUONQ02::MUONQ02" ) { h.MU = (*lvs)[i]; }
37  else if ( (*lvs)[i]->GetName() == "TTR_BARREL::TTR_BARREL" ) { h.TTR = (*lvs)[i]; }
38  }
39 
40  const auto& worldVolume = G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume()->GetLogicalVolume();
41  ISFG4GeoHelper::checkVolumeDepth(worldVolume, truthVolLevel);
42  return h;
43  }();
44 
45  AtlasDetDescr::AtlasRegion nextGeoID = truthVolLevel > 1 ? AtlasDetDescr::fAtlasCavern
47 
48  const G4StepPoint *postStep = aStep->GetPostStepPoint();
49  bool leavingG4World = postStep->GetStepStatus()==fWorldBoundary;
50 
51  if ( leavingG4World ) {
52  nextGeoID = AtlasDetDescr::fAtlasCavern;
53  return nextGeoID;
54  }
55 
56  // If in mother volume, use the ISF_GeoIDSvc to resolve the geoID
57  if (G4StepHelper::postStepBranchDepth(aStep)<truthVolLevel) {
58  nextGeoID = getNextGeoIDFromSvc(*postStep, *geoIDSvc);
59  return nextGeoID;
60  }
61 
62  // Ordering inside out (most truth in the ID anyway...)
63  if ( vh.ID==G4StepHelper::getPostStepLogicalVolume(aStep,truthVolLevel) ) {
64  nextGeoID = AtlasDetDescr::fAtlasID;
65  }
66  else if ( vh.CALO==G4StepHelper::getPostStepLogicalVolume(aStep,truthVolLevel) ) {
67  nextGeoID = AtlasDetDescr::fAtlasCalo;
68  }
69  else if ( vh.MU==G4StepHelper::getPostStepLogicalVolume(aStep,truthVolLevel) ) {
70  nextGeoID = AtlasDetDescr::fAtlasMS;
71  }
72  else if ( vh.BP==G4StepHelper::getPostStepLogicalVolume(aStep,truthVolLevel) ) {
73  nextGeoID = (G4StepHelper::postStepBranchDepth(aStep)>truthVolLevel && G4StepHelper::getPostStepLogicalVolumeName(aStep,truthVolLevel+1)=="BeamPipe::BeamPipeCentral")?AtlasDetDescr::fAtlasID:AtlasDetDescr::fAtlasForward;
74  }
75  else if ( vh.TTR==G4StepHelper::getPostStepLogicalVolume(aStep,truthVolLevel) ) {
76  nextGeoID = AtlasDetDescr::fAtlasCavern;
77  }
78  else if (truthVolLevel>0 && G4StepHelper::getPostStepLogicalVolumeName(aStep,truthVolLevel-1).find("CavernInfra")!=std::string::npos) {
79  nextGeoID = AtlasDetDescr::fAtlasCavern;
80  }
81  else {
82  nextGeoID = getNextGeoIDFromSvc(*postStep, *geoIDSvc);
83  }
84 
85  return nextGeoID;
86 }

The documentation for this class was generated from the following files:
egammaParameters::depth
@ depth
pointing depth of the shower as calculated in egammaqgcld
Definition: egammaParamDefs.h:276
AtlasDetDescr::fAtlasForward
@ fAtlasForward
Definition: AtlasRegion.h:34
ID
std::vector< Identifier > ID
Definition: CalibHitIDCheck.h:24
AtlasDetDescr::AtlasRegion
AtlasRegion
Definition: AtlasRegion.h:27
AtlasDetDescr::fUndefinedAtlasRegion
@ fUndefinedAtlasRegion
Definition: AtlasRegion.h:29
G4StepHelper::postStepBranchDepth
int postStepBranchDepth(const G4Step *theStep)
TODO.
Definition: StepHelper.cxx:82
iGeant4::ISFG4GeoHelper::checkVolumeDepth
static bool checkVolumeDepth(G4LogicalVolume *logicalVol, int volLevel, int depth=0)
Definition: ISFG4GeoHelper.cxx:109
AtlasDetDescr::fAtlasMS
@ fAtlasMS
Definition: AtlasRegion.h:36
lumiFormat.i
int i
Definition: lumiFormat.py:92
extractSporadic.h
list h
Definition: extractSporadic.py:97
egammaPIDObs::CALO
const unsigned int CALO
all cuts in calorimeter (including isolation)
Definition: egammaPIDdefsObs.h:919
G4StepHelper::getPostStepLogicalVolume
G4LogicalVolume * getPostStepLogicalVolume(const G4Step *theStep, int iLevel=0)
TODO.
Definition: StepHelper.cxx:58
AtlasDetDescr::fAtlasCavern
@ fAtlasCavern
Definition: AtlasRegion.h:37
AtlasDetDescr::fAtlasID
@ fAtlasID
Definition: AtlasRegion.h:33
AtlasDetDescr::fAtlasCalo
@ fAtlasCalo
Definition: AtlasRegion.h:35
h
ISF::IGeoIDSvc::identifyGeoID
virtual AtlasDetDescr::AtlasRegion identifyGeoID(const Amg::Vector3D &pos) const =0
A static filter that returns the AtlasRegion of the given ISFParticle (position) -> returns ISF::fUnd...
python.PileUpEventType.Cavern
int Cavern
Definition: PileUpEventType.py:5
iGeant4::ISFG4GeoHelper::getNextGeoIDFromSvc
static AtlasDetDescr::AtlasRegion getNextGeoIDFromSvc(const G4StepPoint &postStep, const ISF::IGeoIDSvc &geoIDSvc)
get the next GeoID using only the geoIDSvc
Definition: ISFG4GeoHelper.cxx:88
G4StepHelper::getPostStepLogicalVolumeName
std::string getPostStepLogicalVolumeName(const G4Step *theStep, int iLevel=0)
TODO.
Definition: StepHelper.cxx:62
description
std::string description
glabal timer - how long have I taken so far?
Definition: hcg.cxx:88