BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSSDVolumeExit.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSDebug.hh"
20#include "BDSGlobalConstants.hh"
21#include "BDSHitEnergyDepositionGlobal.hh"
22#include "BDSSDVolumeExit.hh"
23
24#include "globals.hh"
25#include "G4SDManager.hh"
26#include "G4StepPoint.hh"
27#include "G4StepStatus.hh"
28#include "G4Track.hh"
29
30BDSSDVolumeExit::BDSSDVolumeExit(G4String name,
31 G4bool worldExitIn):
32 G4VSensitiveDetector("volume_exit/" + name),
33 colName(name),
34 HCIDve(-1),
35 hits(nullptr)
36{
37 collectionName.insert(name);
38
39 statusToMatch = worldExitIn ? G4StepStatus::fWorldBoundary : fGeomBoundary;
40}
41
42void BDSSDVolumeExit::Initialize(G4HCofThisEvent* HCE)
43{
45 if (HCIDve < 0)
46 {HCIDve = G4SDManager::GetSDMpointer()->GetCollectionID(hits);}
47 HCE->AddHitsCollection(HCIDve, hits);
48
49#ifdef BDSDEBUG
50 G4cout << __METHOD_NAME__ << "Hits Collection ID: " << HCIDve << G4endl;
51#endif
52}
53
54G4bool BDSSDVolumeExit::ProcessHits(G4Step* aStep,
55 G4TouchableHistory* /*th*/)
56{
57 G4StepPoint* postStepPoint = aStep->GetPostStepPoint();
58
59 if (postStepPoint->GetStepStatus() == statusToMatch)
60 {
61 G4double totalEnergy = postStepPoint->GetTotalEnergy();
62 G4double preStepKineticEnergy = aStep->GetPreStepPoint()->GetKineticEnergy();
63 G4double postStepKineticEnergy = postStepPoint->GetKineticEnergy();
64 G4double stepLength = aStep->GetStepLength();
65 G4ThreeVector position = postStepPoint->GetPosition();
66 G4double T = postStepPoint->GetGlobalTime();
67 G4Track* track = aStep->GetTrack();
68 G4int pdgID = track->GetDefinition()->GetPDGEncoding();
69 G4int trackID = track->GetTrackID();
70 G4int parentID = track->GetParentID();
71 G4double weight = track->GetWeight();
72 G4int turnsTaken = BDSGlobalConstants::Instance()->TurnsTaken();
73
75 preStepKineticEnergy,
76 postStepKineticEnergy,
77 stepLength,
78 position.x(),
79 position.y(),
80 position.z(),
81 T,
82 pdgID,
83 trackID,
84 parentID,
85 weight,
86 turnsTaken);
87 hits->insert(hit);
88 return true;
89 }
90 else
91 {return false;}
92}
static BDSGlobalConstants * Instance()
Access method.
Information recorded for a step leaving a volume.
G4int HCIDve
Hits collection ID for volume exit.
G4String colName
Name prepared for collection.
G4StepStatus statusToMatch
World or volume exit status cache.
BDSHitsCollectionEnergyDepositionGlobal * hits
Hits collection.
const char * GetName(int)
Name of element.
Definition: python.cc:38