BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSSDTerminator.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 "BDSSDTerminator.hh"
22
23#include "G4ios.hh"
24#include "G4TouchableHistory.hh"
25#include "G4VTouchable.hh"
26#include "G4Step.hh"
27
28#include <iomanip>
29
31
32
33BDSSDTerminator::BDSSDTerminator(G4String name)
34 :G4VSensitiveDetector(name)
35{
36 moduloEvents = BDSGlobalConstants::Instance()->PrintModuloEvents();
37 moduloTurns = BDSGlobalConstants::Instance()->PrintModuloTurns();
38}
39
40BDSSDTerminator::~BDSSDTerminator()
41{
42 eventNumber = 0;
43}
44
45void BDSSDTerminator::Initialize(G4HCofThisEvent* /*HCE*/)
46{
48 //we don't actually use HCE here as we don't need to log any of the particle info
49}
50
51G4bool BDSSDTerminator::ProcessHits(G4Step* aStep, G4TouchableHistory*)
52{
53 G4Track* theTrack = aStep->GetTrack();
54 G4int parentID = theTrack->GetParentID();
55 G4double trackLength = theTrack->GetTrackLength();
56#ifdef BDSDEBUG
57 G4cout << __METHOD_NAME__ << "parentID: " << parentID << G4endl;
58 G4cout << __METHOD_NAME__ << "track lenth (mm): " << trackLength << G4endl;
59#endif
60
61 if ((parentID == 0) && (trackLength > 1000) && aStep->GetPreStepPoint()->GetStepStatus() == fGeomBoundary)
62 {
63 // parentID == 0 -> primary particle - should only increment turn number for primaries
64 // trackLength > 1000 (mm) -> not due to initial coordinate offset (at least 1 turn)
65 // pre step point status is on a geometry boundary -> particle just entering volume
66 // and not starting in middle
67#ifdef BDSDEBUG
68 G4cout << __METHOD_NAME__ << "Incrementing turn number " << G4endl;
69 G4cout << __METHOD_NAME__ << "Primary particle - incrementing turn number" << G4endl;
70 G4cout << __METHOD_NAME__ << "Track length is : " << trackLength << " m" << G4endl;
71 G4cout << __METHOD_NAME__ << "Turn number is : " << BDSGlobalConstants::Instance()->TurnsTaken() << G4endl;
72#endif
73 if (eventNumber % moduloEvents == 0)
74 {// feedback on this event
75 G4int turnstaken = BDSGlobalConstants::Instance()->TurnsTaken();
76 if (turnstaken % moduloTurns == 0)
77 {
78 std::ios_base::fmtflags ff = G4cout.flags();// save flags since they are changed
79 G4cout << "Turn: " << std::right << std::setw(4) << std::fixed
80 << turnstaken << " / " << std::left
81 << BDSGlobalConstants::Instance()->TurnsToTake() << G4endl;
82
83 G4cout.flags(ff);// reset flags
84 }
85 }
86
88
89#ifdef BDSDEBUG
90 G4cout << __METHOD_NAME__ << "New turn number : " << BDSGlobalConstants::Instance()->TurnsTaken() << G4endl;
91#endif
92 }
93#ifdef BDSDEBUG
94 else
95 {G4cout << __METHOD_NAME__ << "not incrementing turn number" << G4endl;}
96#endif
97 return true;
98}
99
100
101void BDSSDTerminator::EndOfEvent(G4HCofThisEvent* /*HCE*/)
102{
104}
void IncrementTurnNumber()
Setter.
static BDSGlobalConstants * Instance()
Access method.
void ResetTurnNumber()
Setter.
G4int moduloTurns
Cache of print turn number on these turns.
static G4int eventNumber
Externally accessible counter for event number. Set in BeginOfEventAction.
G4int moduloEvents
Cache of print turn number on these events.