Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HEPEvtInterface.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4HEPEvtInterface class implementation
27//
28// Author: Makoto Asai, 1997
29// --------------------------------------------------------------------
30
31#include "G4HEPEvtInterface.hh"
32
33#include "G4Types.hh"
34#include "G4SystemOfUnits.hh"
35
36#include "G4ios.hh"
37#include "G4PrimaryVertex.hh"
38#include "G4PrimaryParticle.hh"
39#include "G4HEPEvtParticle.hh"
40#include "G4Event.hh"
41
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}
61
63{
64 G4int NHEP = 0; // number of entries
65 if (inputFile.is_open())
66 {
67 inputFile >> NHEP;
68 }
69 else
70 {
71 G4Exception("G4HEPEvtInterface::G4HEPEvtInterface","Event0201",
72 FatalException, "G4HEPEvtInterface:: cannot open file.");
73 }
74 if( inputFile.eof() )
75 {
76 G4Exception("G4HEPEvtInterface::GeneratePrimaryVertex", "Event0202",
78 "End-Of-File: HEPEvt input file -- no more event to read!");
79 return;
80 }
81
82 if(vLevel > 0)
83 {
84 G4cout << "G4HEPEvtInterface - reading " << NHEP
85 << " HEPEvt particles from " << fileName << "." << G4endl;
86 }
87 for( G4int IHEP=0; IHEP<NHEP; ++IHEP )
88 {
89 G4int ISTHEP; // status code
90 G4int IDHEP; // PDG code
91 G4int JDAHEP1; // first daughter
92 G4int JDAHEP2; // last daughter
93 G4double PHEP1; // px in GeV
94 G4double PHEP2; // py in GeV
95 G4double PHEP3; // pz in GeV
96 G4double PHEP5; // mass in GeV
97
98 inputFile >> ISTHEP >> IDHEP >> JDAHEP1 >> JDAHEP2
99 >> PHEP1 >> PHEP2 >> PHEP3 >> PHEP5;
100 if( inputFile.eof() )
101 {
102 G4Exception("G4HEPEvtInterface::GeneratePrimaryVertex", "Event0203",
104 "Unexpected End-Of-File in the middle of an event");
105 }
106 if(vLevel > 1)
107 {
108 G4cout << " " << ISTHEP << " " << IDHEP << " " << JDAHEP1
109 << " " << JDAHEP2 << " " << PHEP1 << " " << PHEP2
110 << " " << PHEP3 << " " << PHEP5 << G4endl;
111 }
112
113 // Create G4PrimaryParticle object
114 //
115 auto* particle = new G4PrimaryParticle( IDHEP );
116 particle->SetMass( PHEP5*GeV );
117 particle->SetMomentum(PHEP1*GeV, PHEP2*GeV, PHEP3*GeV );
118
119 // Create G4HEPEvtParticle object
120 //
121 auto* hepParticle
122 = new G4HEPEvtParticle( particle, ISTHEP, JDAHEP1, JDAHEP2 );
123
124 // Store
125 //
126 HPlist.push_back( hepParticle );
127 }
128
129 // Check if there is at least one particle
130 //
131 if( HPlist.empty() ) return;
132
133 // Make connection between daughter particles decayed from the same mother
134 //
135 for(auto & i : HPlist)
136 {
137 if( i->GetJDAHEP1() > 0 ) // it has daughters
138 {
139 G4int jda1 = i->GetJDAHEP1()-1; // FORTRAN index starts from 1
140 G4int jda2 = i->GetJDAHEP2()-1; // but C++ starts from 0.
141 G4PrimaryParticle* mother = i->GetTheParticle();
142 for( G4int j=jda1; j<=jda2; ++j )
143 {
144 G4PrimaryParticle* daughter = HPlist[j]->GetTheParticle();
145 if(HPlist[j]->GetISTHEP()>0)
146 {
147 mother->SetDaughter( daughter );
148 HPlist[j]->Done();
149 }
150 }
151 }
152 }
153
154 // Create G4PrimaryVertex object
155 //
157
158 // Put initial particles to the vertex
159 //
160 for(auto & ii : HPlist)
161 {
162 if( ii->GetISTHEP() > 0 ) // ISTHEP of daughters had been
163 // set to negative
164 {
165 G4PrimaryParticle* initialParticle = ii->GetTheParticle();
166 vertex->SetPrimary( initialParticle );
167 }
168 }
169
170 // Clear G4HEPEvtParticles
171 //
172 for(auto & iii : HPlist)
173 {
174 delete iii;
175 }
176 HPlist.clear();
177
178 // Put the vertex to G4Event object
179 //
180 evt->AddPrimaryVertex( vertex );
181}
@ FatalException
@ RunMustBeAborted
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition G4Event.hh:126
G4HEPEvtInterface(const char *evfile, G4int vl=0)
void GeneratePrimaryVertex(G4Event *evt) override
void SetDaughter(G4PrimaryParticle *np)
G4ThreeVector particle_position