BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBunchFileBased.hh
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#ifndef BDSBUNCHFILEBASED_H
20#define BDSBUNCHFILEBASED_H
21
22#include "BDSBunch.hh"
23
24#include "G4String.hh"
25
38{
39public:
40 explicit BDSBunchFileBased(const G4String& distributionName);
41 virtual ~BDSBunchFileBased();
46
49 virtual void BeginOfRunAction(G4int numberOfEvents,
50 G4bool batchMode);
51
53 unsigned long long int NOriginalEvents() const {return nOriginalEvents;}
54 unsigned long long int NEventsInFile() const {return nEventsInFile;}
55 unsigned long long int NEventsInFileSkipped() const {return nEventsInFileSkipped;}
56 G4int DistrFileLoopNTimes() const {return distrFileLoopNTimes;}
58
59 void SetNEventsInFile(unsigned long long int nEventsInFileIn) {nEventsInFile = nEventsInFileIn;}
60 void SetNOriginalEvents(unsigned long long int nOriginalEventsIn) {nOriginalEvents = nOriginalEventsIn;}
61
62 void IncrementNEventsInFileSkipped() {nEventsInFileSkipped += 1;}
63 void IncrementNEventsInFileSkipped(unsigned long long int plus) {nEventsInFileSkipped += plus;}
64
66 virtual void SetOptions(const BDSParticleDefinition* beamParticle,
67 const GMAD::Beam& beam,
68 const BDSBunchType& distrType,
69 G4Transform3D beamlineTransformIn = G4Transform3D::Identity,
70 const G4double beamlineS = 0);
71
72protected:
73 unsigned long long int nOriginalEvents;
74 unsigned long long int nEventsInFile;
75 unsigned long long int nEventsInFileSkipped;
76 G4bool distrFileLoop;
77 G4int distrFileLoopNTimes;
78};
79
80#endif
An intermediate layer for any bunch classes that are file based.
BDSBunchFileBased(BDSBunchFileBased &)=delete
Assignment and copy constructor not implemented nor used.
unsigned long long int nEventsInFile
The number of entries in the file loaded.
unsigned long long int NOriginalEvents() const
Accessor.
virtual void BeginOfRunAction(G4int numberOfEvents, G4bool batchMode)
unsigned long long int NEventsInFileSkipped() const
Accessor.
unsigned long long int nEventsInFileSkipped
Number that are skipped as we go through the file due to filters.
G4int DistrFileLoopNTimes() const
Accessor.
unsigned long long int nOriginalEvents
nOriginalEvents from upstream file if skimmed - need to pass through.
virtual void SetOptions(const BDSParticleDefinition *beamParticle, const GMAD::Beam &beam, const BDSBunchType &distrType, G4Transform3D beamlineTransformIn=G4Transform3D::Identity, const G4double beamlineS=0)
Pull out the relevant options and then pass through to BDSBunch::SetOptions().
BDSBunchFileBased & operator=(const BDSBunchFileBased &)=delete
Assignment and copy constructor not implemented nor used.
unsigned long long int NEventsInFile() const
Accessor.
The base class for bunch distribution generators.
Definition: BDSBunch.hh:47
G4double beamlineS
Beamline initial S position.
Definition: BDSBunch.hh:219
Wrapper for particle definition.
Improve type-safety of native enum data type in C++.
Beam class.
Definition: beam.h:44