BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSPrimaryGeneratorFileSampler.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 BDSPRIMARYGENERATORFILESAMPLER_H
20#define BDSPRIMARYGENERATORFILESAMPLER_H
21
22#include "BDSPrimaryGeneratorFile.hh"
23
24#include "G4RotationMatrix.hh"
25#include "G4String.hh"
26#include "G4ThreeVector.hh"
27#include "G4Types.hh"
28
29#include <vector>
30
33class G4Event;
34class G4PrimaryParticle;
35class G4PrimaryVertex;
36class G4VSolid;
37
45{
46public:
53 BDSPrimaryGeneratorFileSampler(const G4String& distrType,
54 const G4String& fileNameIn,
56 G4bool loopFileIn,
57 G4bool removeUnstableWithoutDecayIn = true,
58 G4bool warnAboutSkippedParticlesIn = true);
60
62 virtual void GeneratePrimaryVertex(G4Event* anEvent);
63
65 virtual void RecreateAdvanceToEvent(G4int eventOffset);
66
68 {
69 G4ThreeVector xyz;
70 G4PrimaryParticle* vertex;
71 };
72
75 void SkipEvents(G4long nEventsToSkip);
76
77protected:
79 void ReadSingleEvent(G4long index, G4Event* anEvent);
80
82 //void HepMC2G4(const HepMC3::GenEvent* hepmcevt, G4Event* g4event);
83
84 void ReadPrimaryParticlesFloat(G4long index);
85 void ReadPrimaryParticlesDouble(G4long index);
86
87private:
89 G4String fileName;
90 G4String samplerName;
91 G4bool removeUnstableWithoutDecay;
92 G4bool warnAboutSkippedParticles;
93 G4RotationMatrix referenceBeamMomentumOffset;
94
96 std::vector<DisplacedVertex> vertices;
97
98 std::vector<G4PrimaryVertex*> currentVertices;
99};
100
101#endif
A wrapper of BDSBunch to include a filter for the events loaded by an event generator.
Loader of ROOT Event output for receating events.
Loader to read a specific sampler from a BDSIM ROOT output file.
virtual void RecreateAdvanceToEvent(G4int eventOffset)
Advance to the correct event number in the file for recreation.
void ReadSingleEvent(G4long index, G4Event *anEvent)
Read sampler hits and put into primary vertices if they pass filters.
std::vector< DisplacedVertex > vertices
Used for transiently loading information.
BDSPrimaryGeneratorFileSampler()=delete
Do not require default constructor.
void ReadPrimaryParticlesFloat(G4long index)
Conversion from HepMC::GenEvent to G4Event.
virtual void GeneratePrimaryVertex(G4Event *anEvent)
Read the next non-empty sampler entry from the file.
Common interface for any primary generators that are file based.