BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSTerminatorUserLimits.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
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 "BDSTerminatorUserLimits.hh"
22
23#include "globals.hh" // geant4 types / globals
24#include "G4UserLimits.hh"
25
26#include <limits>
27
29 G4double utrakMax,
30 G4double utimeMax,
31 G4double uekinMin,
32 G4double urangMin):
33 G4UserLimits(ustepMax,
34 utrakMax,
35 utimeMax,
36 uekinMin,
37 urangMin),
38 keeprunningEK(0.0),
39 stoprunningEK(std::numeric_limits<double>::max()),
40 turnsToTake(BDSGlobalConstants::Instance()->TurnsToTake())
41{}
42
44 G4double ustepMax,
45 G4double utrakMax,
46 G4double utimeMax,
47 G4double uekinMin,
48 G4double urangMin):
49 G4UserLimits(type,
50 ustepMax,
51 utrakMax,
52 utimeMax,
53 uekinMin,
54 urangMin),
55 keeprunningEK(0.0),
56 stoprunningEK(std::numeric_limits<double>::max()),
57 turnsToTake(BDSGlobalConstants::Instance()->TurnsToTake())
58{}
59
60inline G4double BDSTerminatorUserLimits::GetUserMinEkine(const G4Track& /*trk*/)
61{
62 // does the number of turns passed == number of turns to take
63 G4int turnsTaken = BDSGlobalConstants::Instance()->TurnsTaken();
64#ifdef BDSDEBUG
65 G4cout << __METHOD_NAME__ << "turns taken : " << turnsTaken << G4endl;
66#endif
67 if (turnsTaken >= turnsToTake)
68 {
69#ifdef BDSDEBUG
70 G4cout << __METHOD_NAME__ << "Requested number of turns completed - stopping all particles" << G4endl;
71#endif
72 return stoprunningEK;
73 } // yes - stop: return DBL_MAX eV so no particles will be tracked
74 else
75 {
76 return keeprunningEK;
77 } // no - keep going: return 0 eV so all are tracked
78}
79
80
A class that holds global options and constants.
static BDSGlobalConstants * Instance()
Access method.
virtual G4double GetUserMinEkine(const G4Track &)
G4double keeprunningEK
Minimum energy particle must have to keep going.
G4double stoprunningEK
Same, so everything < DBL_MAX so everything stopped.
BDSTerminatorUserLimits(G4double ustepMax=std::numeric_limits< double >::max(), G4double utrakMax=std::numeric_limits< double >::max(), G4double utimeMax=std::numeric_limits< double >::max(), G4double uekinMin=0., G4double urangMin=0.)