Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HEPEvtInterface Class Reference

#include <G4HEPEvtInterface.hh>

+ Inheritance diagram for G4HEPEvtInterface:

Public Member Functions

 G4HEPEvtInterface (const char *evfile, G4int vl=0)
 
 ~G4HEPEvtInterface ()
 
void GeneratePrimaryVertex (G4Event *evt)
 
- Public Member Functions inherited from G4VPrimaryGenerator
 G4VPrimaryGenerator ()
 
virtual ~G4VPrimaryGenerator ()
 
virtual void GeneratePrimaryVertex (G4Event *evt)=0
 
G4ThreeVector GetParticlePosition ()
 
G4double GetParticleTime ()
 
void SetParticlePosition (G4ThreeVector aPosition)
 
void SetParticleTime (G4double aTime)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPrimaryGenerator
static G4bool CheckVertexInsideWorld (const G4ThreeVector &pos)
 
- Protected Attributes inherited from G4VPrimaryGenerator
G4ThreeVector particle_position
 
G4double particle_time = 0.0
 

Detailed Description

Definition at line 78 of file G4HEPEvtInterface.hh.

Constructor & Destructor Documentation

◆ G4HEPEvtInterface()

G4HEPEvtInterface::G4HEPEvtInterface ( const char *  evfile,
G4int  vl = 0 
)

Definition at line 42 of file G4HEPEvtInterface.cc.

43 : vLevel(vl)
44{
45 inputFile.open((char*)evfile);
46 if (inputFile.is_open())
47 {
48 fileName = evfile;
49 if(vl>0)
50 G4cout << "G4HEPEvtInterface - " << fileName << " is open." << G4endl;
51 }
52 else
53 {
54 G4Exception("G4HEPEvtInterface::G4HEPEvtInterface","Event0201",
55 FatalException, "G4HEPEvtInterface:: cannot open file.");
56 }
57 G4ThreeVector zero;
58 particle_position = zero;
59 particle_time = 0.0;
60}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
G4ThreeVector particle_position

◆ ~G4HEPEvtInterface()

G4HEPEvtInterface::~G4HEPEvtInterface ( )

Definition at line 62 of file G4HEPEvtInterface.cc.

63{
64}

Member Function Documentation

◆ GeneratePrimaryVertex()

void G4HEPEvtInterface::GeneratePrimaryVertex ( G4Event evt)
virtual

Implements G4VPrimaryGenerator.

Definition at line 66 of file G4HEPEvtInterface.cc.

67{
68 G4int NHEP = 0; // number of entries
69 if (inputFile.is_open())
70 {
71 inputFile >> NHEP;
72 }
73 else
74 {
75 G4Exception("G4HEPEvtInterface::G4HEPEvtInterface","Event0201",
76 FatalException, "G4HEPEvtInterface:: cannot open file.");
77 }
78 if( inputFile.eof() )
79 {
80 G4Exception("G4HEPEvtInterface::GeneratePrimaryVertex", "Event0202",
82 "End-Of-File: HEPEvt input file -- no more event to read!");
83 return;
84 }
85
86 if(vLevel > 0)
87 {
88 G4cout << "G4HEPEvtInterface - reading " << NHEP
89 << " HEPEvt particles from " << fileName << "." << G4endl;
90 }
91 for( G4int IHEP=0; IHEP<NHEP; ++IHEP )
92 {
93 G4int ISTHEP; // status code
94 G4int IDHEP; // PDG code
95 G4int JDAHEP1; // first daughter
96 G4int JDAHEP2; // last daughter
97 G4double PHEP1; // px in GeV
98 G4double PHEP2; // py in GeV
99 G4double PHEP3; // pz in GeV
100 G4double PHEP5; // mass in GeV
101
102 inputFile >> ISTHEP >> IDHEP >> JDAHEP1 >> JDAHEP2
103 >> PHEP1 >> PHEP2 >> PHEP3 >> PHEP5;
104 if( inputFile.eof() )
105 {
106 G4Exception("G4HEPEvtInterface::GeneratePrimaryVertex", "Event0203",
108 "Unexpected End-Of-File in the middle of an event");
109 }
110 if(vLevel > 1)
111 {
112 G4cout << " " << ISTHEP << " " << IDHEP << " " << JDAHEP1
113 << " " << JDAHEP2 << " " << PHEP1 << " " << PHEP2
114 << " " << PHEP3 << " " << PHEP5 << G4endl;
115 }
116
117 // Create G4PrimaryParticle object
118 //
119 G4PrimaryParticle* particle = new G4PrimaryParticle( IDHEP );
120 particle->SetMass( PHEP5*GeV );
121 particle->SetMomentum(PHEP1*GeV, PHEP2*GeV, PHEP3*GeV );
122
123 // Create G4HEPEvtParticle object
124 //
125 G4HEPEvtParticle* hepParticle
126 = new G4HEPEvtParticle( particle, ISTHEP, JDAHEP1, JDAHEP2 );
127
128 // Store
129 //
130 HPlist.push_back( hepParticle );
131 }
132
133 // Check if there is at least one particle
134 //
135 if( HPlist.size() == 0 ) return;
136
137 // Make connection between daughter particles decayed from the same mother
138 //
139 for( std::size_t i=0; i<HPlist.size(); ++i )
140 {
141 if( HPlist[i]->GetJDAHEP1() > 0 ) // it has daughters
142 {
143 G4int jda1 = HPlist[i]->GetJDAHEP1()-1; // FORTRAN index starts from 1
144 G4int jda2 = HPlist[i]->GetJDAHEP2()-1; // but C++ starts from 0.
145 G4PrimaryParticle* mother = HPlist[i]->GetTheParticle();
146 for( G4int j=jda1; j<=jda2; ++j )
147 {
148 G4PrimaryParticle* daughter = HPlist[j]->GetTheParticle();
149 if(HPlist[j]->GetISTHEP()>0)
150 {
151 mother->SetDaughter( daughter );
152 HPlist[j]->Done();
153 }
154 }
155 }
156 }
157
158 // Create G4PrimaryVertex object
159 //
161
162 // Put initial particles to the vertex
163 //
164 for( std::size_t ii=0; ii<HPlist.size(); ++ii )
165 {
166 if( HPlist[ii]->GetISTHEP() > 0 ) // ISTHEP of daughters had been
167 // set to negative
168 {
169 G4PrimaryParticle* initialParticle = HPlist[ii]->GetTheParticle();
170 vertex->SetPrimary( initialParticle );
171 }
172 }
173
174 // Clear G4HEPEvtParticles
175 //
176 for(std::size_t iii=0; iii<HPlist.size(); ++iii)
177 {
178 delete HPlist[iii];
179 }
180 HPlist.clear();
181
182 // Put the vertex to G4Event object
183 //
184 evt->AddPrimaryVertex( vertex );
185}
@ RunMustBeAborted
double G4double
Definition: G4Types.hh:83
int G4int
Definition: G4Types.hh:85
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:121
void SetMomentum(G4double px, G4double py, G4double pz)
void SetMass(G4double mas)
void SetDaughter(G4PrimaryParticle *np)
void SetPrimary(G4PrimaryParticle *pp)

The documentation for this class was generated from the following files: