BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSHepMC3Reader.hh
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#ifdef USE_HEPMC3
20
21#ifndef BDSHEPMC3READER_H
22#define BDSHEPMC3READER_H
23
24#include "BDSEventGeneratorFileType.hh"
25
26#include "globals.hh"
27#include "G4RotationMatrix.hh"
28#include "G4ThreeVector.hh"
29#include "G4VPrimaryGenerator.hh"
30
32class G4Event;
33class G4VSolid;
34
35namespace HepMC3
36{
37 class GenEvent;
38 class Reader;
39}
40
55class BDSHepMC3Reader: public G4VPrimaryGenerator
56{
57public:
59 BDSHepMC3Reader() = delete;
64 BDSHepMC3Reader(const G4String& distrType,
65 const G4String& fileNameIn,
67 G4bool removeUnstableWithoutDecayIn = true,
68 G4bool warnAboutSkippedParticlesIn = true);
69 virtual ~BDSHepMC3Reader();
70
72 inline HepMC3::GenEvent* GetHepMCGenEvent() const {return hepmcEvent;}
73
74 // The default behavior is that a single HepMC event generated by
75 // GenerateHepMCEvent() will be converted to G4Event through HepMC2G4().
76 virtual void GeneratePrimaryVertex(G4Event* anEvent);
77
79 virtual void RecreateAdvanceToEvent(G4int eventOffset);
80
81protected:
83 void OpenFile();
84
87 void CloseFile();
88
90 void ReadSingleEvent();
91
93 void HepMC2G4(const HepMC3::GenEvent* hepmcevt, G4Event* g4event);
94
95 // We have to take care for the position of primaries because
96 // primary vertices outside the world volume would give rise to a G4Exception.
97 virtual G4bool VertexInsideWorld(const G4ThreeVector& pos) const;
98
99 // Note that the life of HepMC event object must be handled by users.
100 // In the default implementation, a current HepMC event will be
101 // deleted at GeneratePrimaryVertex() in the next event.
102 HepMC3::GenEvent* hepmcEvent; // (care for single event case only)
103
104private:
105 HepMC3::Reader* reader;
106 G4String fileName;
108 G4bool removeUnstableWithoutDecay;
109 G4bool warnAboutSkippedParticles;
111 G4RotationMatrix referenceBeamMomentumOffset;
112 mutable G4VSolid* worldSolid;
113};
114
115#endif
116#endif
A wrapper of BDSBunch to include a filter for the events loaded by an event generator.
Loader to use any HepMC3 compatible file.
BDSHepMC3Reader()=delete
Do not require default constructor.
virtual void RecreateAdvanceToEvent(G4int eventOffset)
Advance to the correct event number in the file for recreation.
void OpenFile()
Construct the member "reader" and open the file for reading.
void ReadSingleEvent()
Clear the hepmcEvent object, reallocate and read a single event and fill that member.
void HepMC2G4(const HepMC3::GenEvent *hepmcevt, G4Event *g4event)
Conversion from HepMC::GenEvent to G4Event.
HepMC3::GenEvent * GetHepMCGenEvent() const
Accessor.