ATLAS Offline Software
Loading...
Searching...
No Matches
ActsGeantFollower.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4
6// GeantFollower.cxx, (c) ATLAS Detector software
8
9#include "ActsGeantFollower.h"
12
13#include "G4Event.hh"
14#include "G4Step.hh"
15#include "G4Material.hh"
16#include "G4StepPoint.hh"
17#include "G4TouchableHistory.hh"
18#include "G4LogicalVolume.hh"
19#include "G4DynamicParticle.hh"
20#include "G4Track.hh"
21#include "G4VSensitiveDetector.hh"
22
24
26{
27 m_helper->beginEvent();
28}
29
31{
32 m_helper->endEvent();
33}
34
36{
37 if(m_helper.retrieve()!=StatusCode::SUCCESS)
38 {
39 G4ExceptionDescription description;
40 description << "Cannot retrieve ActsGeantFollower helper";
41 G4Exception("ActsGeantFollower", "ActsGeantFollower1", FatalException, description);
42 return;
43 }
44 return;
45}
46
48{
49 // kill secondaries and low momentum particles
50 if (aStep->GetTrack()->GetParentID() || aStep->GetPreStepPoint()->GetMomentum().mag()<500 )
51 {
52 std::cout << "low pt" << std::endl;
53 aStep->GetTrack()->SetTrackStatus(fStopAndKill);
54 return;
55 }
56
57 // kill Particles outiside the tracking volume
58 if (aStep->GetPreStepPoint()->GetPosition().z()>3000 || sqrt(aStep->GetPreStepPoint()->GetPosition().x()*aStep->GetPreStepPoint()->GetPosition().x()+aStep->GetPreStepPoint()->GetPosition().y()*aStep->GetPreStepPoint()->GetPosition().y())>1050 )
59 {
60 std::cout << "out" << std::endl;
61 aStep->GetTrack()->SetTrackStatus(fStopAndKill);
62 return;
63 }
64
65 // get the prestep point and follow this guy
66 G4StepPoint * g4PreStep = aStep->GetPreStepPoint();
67 G4ThreeVector g4Momentum = g4PreStep->GetMomentum();
68 G4ThreeVector g4Position = g4PreStep->GetPosition();
69
70 G4Track* g4Track = aStep->GetTrack();
71 const G4DynamicParticle* g4DynParticle = g4Track->GetDynamicParticle();
72
73 // the material information
74 const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
75 if(ATH_LIKELY(touchHist))
76 {
77 // G4LogicalVolume
78 const G4LogicalVolume *lv= touchHist->GetVolume()->GetLogicalVolume();
79 if(ATH_LIKELY(lv))
80 {
81 const G4Material *mat = lv->GetMaterial();
82 // the step information
83 double steplength = aStep->GetStepLength();
84 // the position information
85 double X0 = mat->GetRadlen();
86 // update the track follower when a sensor is encountered
87 // bool isSensitive = (lv->GetSensitiveDetector() != nullptr);
88 bool isSensitive = true;
89 m_helper->trackParticle(g4Position, g4Momentum, g4DynParticle->GetPDGcode(), g4DynParticle->GetCharge(), steplength, X0, isSensitive);
90
91 }
92 else
93 {
94 G4ExceptionDescription description;
95 description << "ActsGeantFollower::SteppingAction NULL G4LogicalVolume pointer.";
96 G4Exception("ActsGeantFollower", "ActsGeantFollower2", FatalException, description);
97 }
98 }
99 else
100 {
101 G4ExceptionDescription description;
102 description << "ActsGeantFollower::SteppingAction NULL G4TouchableHistory pointer.";
103 G4Exception("ActsGeantFollower", "ActsGeantFollower3", FatalException, description);
104 }
105 return;
106}
#define ATH_LIKELY(x)
virtual void EndOfEventAction(const G4Event *) override
ToolHandle< IActsGeantFollowerHelper > m_helper
virtual void BeginOfEventAction(const G4Event *) override
virtual void BeginOfRunAction(const G4Run *) override
virtual void UserSteppingAction(const G4Step *) override
std::string description
glabal timer - how long have I taken so far?
Definition hcg.cxx:91