ATLAS Offline Software
Loading...
Searching...
No Matches
TestBoundariesUserAction.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
7#include "G4ParticleDefinition.hh"
8#include "G4DynamicParticle.hh"
9#include "G4TouchableHistory.hh"
10#include "G4Step.hh"
11
12#include "TFile.h"
13#include "TTree.h"
14
15#include "GaudiKernel/Bootstrap.h"
16#include "GaudiKernel/ISvcLocator.h"
17#include "GaudiKernel/IMessageSvc.h"
18
19namespace G4UA
20{
21
23 : AthMessaging(Gaudi::svcLocator()->service< IMessageSvc >( "MessageSvc" ),
24 "TestBoundariesUserAction")
25 {
26 }
27
29 {
30 m_file = TFile::Open("points.root","RECREATE");
31 ATH_MSG_INFO("Open file points.root, create tree");
32 m_file->cd();
33 m_tree = new TTree("points","points");
34 m_tree->Branch("x",&m_data.x,"x/F");
35 m_tree->Branch("y",&m_data.y,"y/F");
36 m_tree->Branch("z",&m_data.z,"z/F");
37 m_tree->Branch("volume",&m_data.volume,"volume/I");
38 m_tree->Branch("enter",&m_data.enter,"enter/O");
39 m_tree->Branch("exit",&m_data.exit,"exit/O");
40 m_tree->Branch("leave",&m_data.leave,"leave/O");
41 }
42
44 {
45 ATH_MSG_INFO("Writing file points.root");
46 m_file->cd();
47 m_tree->Write();
48 m_file->Close();
49 delete m_file;
50 m_file=nullptr;
51 m_tree=nullptr;
52 }
53
55 {
56 G4StepPoint* preStep = aStep->GetPreStepPoint();
57 G4StepPoint* postStep = aStep->GetPostStepPoint();
58
59 if (!preStep || !postStep) {
60 ATH_MSG_ERROR("TestBoundariesUserAction Error something missing");
61 return;
62 }
63
64 const G4TouchableHistory* preTH = dynamic_cast<const G4TouchableHistory*>(preStep->GetTouchable());
65 const G4TouchableHistory* postTH = dynamic_cast<const G4TouchableHistory*>(postStep->GetTouchable());
66 G4LogicalVolume* preLV = preStep->GetPhysicalVolume()->GetLogicalVolume();
67
68 if (preTH) {
69 int preN=preTH->GetHistoryDepth();
70 if (preN>0) preLV = preTH->GetVolume(preN-1)->GetLogicalVolume();
71 }
72
73 G4LogicalVolume* postLV{};
74 if (postStep->GetPhysicalVolume())
75 postLV = postStep->GetPhysicalVolume()->GetLogicalVolume();
76 if (postTH) {
77 int postN = postTH->GetHistoryDepth();
78 if (postN>0) postLV = postTH->GetVolume(postN-1)->GetLogicalVolume();
79 }
80
81 if (postLV!=preLV) {
82 const G4ThreeVector& pos = aStep->GetPostStepPoint()->GetPosition();
83 m_data.Set(pos.x(),pos.y(),pos.z(),0);
84 if (preLV) {
85 const std::string& preStepPointName = preLV->GetName();
86 SMap::const_iterator it = m_sel.find(preStepPointName);
87 if (it!=m_sel.end()) {
88 m_data.volume = it->second;
89 m_data.Reset();
90 m_data.exit = true;
91 if (m_tree) m_tree->Fill();
92 }
93 }
94 if (postLV) {
95 SMap::const_iterator it = m_sel.find( postLV->GetName() );
96 if (it!=m_sel.end()) {
97 m_data.volume = it->second;
98 m_data.Reset();
99 m_data.enter = true;
100 if (m_tree) m_tree->Fill();
101 }
102 }
103 else {
104 m_data.volume = 0;
105 m_data.Reset();
106 m_data.leave = true;
107 if (m_tree) m_tree->Fill();
108 }
109 }
110 }
111
112} // namespace G4UA
#define ATH_MSG_ERROR(x)
#define ATH_MSG_INFO(x)
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
virtual void UserSteppingAction(const G4Step *) override
virtual void EndOfRunAction(const G4Run *) override
virtual void BeginOfRunAction(const G4Run *) override
struct G4UA::TestBoundariesUserAction::Data m_data
=============================================================================