Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
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//
27// $Id$
28//
29//
30// --------------------------------------------------------------------
31
32#include "G4HEPEvtInterface.hh"
33
34#include "G4Types.hh"
35#include "G4SystemOfUnits.hh"
36
37#include "G4ios.hh"
38#include "G4PrimaryVertex.hh"
39#include "G4PrimaryParticle.hh"
40#include "G4HEPEvtParticle.hh"
41#include "G4Event.hh"
42
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}
58
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}
74
76{;}
77
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}
164
@ JustWarning
@ FatalException
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void AddPrimaryVertex(G4PrimaryVertex *aPrimaryVertex)
Definition: G4Event.hh:142
G4HEPEvtInterface(char *evfile)
void GeneratePrimaryVertex(G4Event *evt)
void SetMomentum(G4double px, G4double py, G4double pz)
void SetMass(G4double mas)
void SetDaughter(G4PrimaryParticle *np)
void SetPrimary(G4PrimaryParticle *pp)
const char * data() const
G4ThreeVector particle_position
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41