Geant4 9.6.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 (char *evfile)
 
 G4HEPEvtInterface (G4String evfile)
 
 ~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
 

Detailed Description

Definition at line 76 of file G4HEPEvtInterface.hh.

Constructor & Destructor Documentation

◆ G4HEPEvtInterface() [1/2]

G4HEPEvtInterface::G4HEPEvtInterface ( char *  evfile)

Definition at line 43 of file G4HEPEvtInterface.cc.

44{
45 inputFile.open(evfile);
46 if (inputFile) {
47 fileName = evfile;
48 }
49 else {
50 G4Exception("G4HEPEvtInterface::G4HEPEvtInterface","Event0201",FatalException,
51 "G4HEPEvtInterface:: cannot open file.");
52 }
53 G4ThreeVector zero;
54 particle_position = zero;
55 particle_time = 0.0;
56
57}
@ FatalException
G4ThreeVector particle_position
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

◆ G4HEPEvtInterface() [2/2]

G4HEPEvtInterface::G4HEPEvtInterface ( G4String  evfile)

Definition at line 59 of file G4HEPEvtInterface.cc.

60{
61 const char* fn = evfile.data();
62 inputFile.open((char*)fn);
63 if (inputFile) {
64 fileName = evfile;
65 }
66 else {
67 G4Exception("G4HEPEvtInterface::G4HEPEvtInterface","Event0201",FatalException,
68 "G4HEPEvtInterface:: cannot open file.");
69 }
70 G4ThreeVector zero;
71 particle_position = zero;
72 particle_time = 0.0;
73}
const char * data() const

◆ ~G4HEPEvtInterface()

G4HEPEvtInterface::~G4HEPEvtInterface ( )

Definition at line 75 of file G4HEPEvtInterface.cc.

76{;}

Member Function Documentation

◆ GeneratePrimaryVertex()

void G4HEPEvtInterface::GeneratePrimaryVertex ( G4Event evt)
virtual

Implements G4VPrimaryGenerator.

Definition at line 78 of file G4HEPEvtInterface.cc.

79{
80 G4int NHEP; // number of entries
81 inputFile >> NHEP;
82 if( inputFile.eof() )
83 {
84 G4Exception("G4HEPEvtInterface::GeneratePrimaryVertex","Event0202",
85 JustWarning,"End-Of-File : HEPEvt input file");
86 return;
87 }
88
89 for( G4int IHEP=0; IHEP<NHEP; IHEP++ )
90 {
91 G4int ISTHEP; // status code
92 G4int IDHEP; // PDG code
93 G4int JDAHEP1; // first daughter
94 G4int JDAHEP2; // last daughter
95 G4double PHEP1; // px in GeV
96 G4double PHEP2; // py in GeV
97 G4double PHEP3; // pz in GeV
98 G4double PHEP5; // mass in GeV
99
100 inputFile >> ISTHEP >> IDHEP >> JDAHEP1 >> JDAHEP2
101 >> PHEP1 >> PHEP2 >> PHEP3 >> PHEP5;
102
103 // create G4PrimaryParticle object
104 G4PrimaryParticle* particle
105 = new G4PrimaryParticle( IDHEP );
106 particle->SetMass( PHEP5*GeV );
107 particle->SetMomentum(PHEP1*GeV, PHEP2*GeV, PHEP3*GeV );
108
109 // create G4HEPEvtParticle object
110 G4HEPEvtParticle* hepParticle
111 = new G4HEPEvtParticle( particle, ISTHEP, JDAHEP1, JDAHEP2 );
112
113 // Store
114 HPlist.push_back( hepParticle );
115 }
116
117 // check if there is at least one particle
118 if( HPlist.size() == 0 ) return;
119
120 // make connection between daughter particles decayed from
121 // the same mother
122 for( size_t i=0; i<HPlist.size(); i++ )
123 {
124 if( HPlist[i]->GetJDAHEP1() > 0 ) // it has daughters
125 {
126 G4int jda1 = HPlist[i]->GetJDAHEP1()-1; // FORTRAN index starts from 1
127 G4int jda2 = HPlist[i]->GetJDAHEP2()-1; // but C++ starts from 0.
128 G4PrimaryParticle* mother = HPlist[i]->GetTheParticle();
129 for( G4int j=jda1; j<=jda2; j++ )
130 {
131 G4PrimaryParticle* daughter = HPlist[j]->GetTheParticle();
132 if(HPlist[j]->GetISTHEP()>0)
133 {
134 mother->SetDaughter( daughter );
135 HPlist[j]->Done();
136 }
137 }
138 }
139 }
140
141 // create G4PrimaryVertex object
143
144 // put initial particles to the vertex
145 for( size_t ii=0; ii<HPlist.size(); ii++ )
146 {
147 if( HPlist[ii]->GetISTHEP() > 0 ) // ISTHEP of daughters had been
148 // set to negative
149 {
150 G4PrimaryParticle* initialParticle = HPlist[ii]->GetTheParticle();
151 vertex->SetPrimary( initialParticle );
152 }
153 }
154
155 // clear G4HEPEvtParticles
156 //HPlist.clearAndDestroy();
157 for(size_t iii=0;iii<HPlist.size();iii++)
158 { delete HPlist[iii]; }
159 HPlist.clear();
160
161 // Put the vertex to G4Event object
162 evt->AddPrimaryVertex( vertex );
163}
@ JustWarning
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:142
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: