source: trunk/source/event/src/G4HEPEvtInterface.cc@ 1026

Last change on this file since 1026 was 964, checked in by garnier, 17 years ago

update event

File size: 5.0 KB
Line 
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: G4HEPEvtInterface.cc,v 1.11 2006/06/29 18:09:48 gunter Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
29//
30//
31// --------------------------------------------------------------------
32
33#include "G4Types.hh"
34
35#include "G4ios.hh"
36#include "G4HEPEvtInterface.hh"
37#include "G4PrimaryVertex.hh"
38#include "G4PrimaryParticle.hh"
39#include "G4HEPEvtParticle.hh"
40#include "G4Event.hh"
41
42G4HEPEvtInterface::G4HEPEvtInterface(char* evfile)
43{
44 inputFile.open(evfile);
45 if (inputFile) {
46 fileName = evfile;
47 }
48 else {
49 G4Exception("G4HEPEvtInterface:: cannot open file.");
50 }
51 G4ThreeVector zero;
52 particle_position = zero;
53 particle_time = 0.0;
54
55}
56
57G4HEPEvtInterface::G4HEPEvtInterface(G4String evfile)
58{
59 const char* fn = evfile.data();
60 inputFile.open((char*)fn);
61 if (inputFile) {
62 fileName = evfile;
63 }
64 else {
65 G4Exception("G4HEPEvtInterface:: cannot open file.");
66 }
67 G4ThreeVector zero;
68 particle_position = zero;
69 particle_time = 0.0;
70}
71
72G4HEPEvtInterface::~G4HEPEvtInterface()
73{;}
74
75void G4HEPEvtInterface::GeneratePrimaryVertex(G4Event* evt)
76{
77 G4int NHEP; // number of entries
78 inputFile >> NHEP;
79 if( inputFile.eof() )
80 {
81 G4Exception("End-Of-File : HEPEvt input file");
82 return;
83 }
84
85 for( G4int IHEP=0; IHEP<NHEP; IHEP++ )
86 {
87 G4int ISTHEP; // status code
88 G4int IDHEP; // PDG code
89 G4int JDAHEP1; // first daughter
90 G4int JDAHEP2; // last daughter
91 G4double PHEP1; // px in GeV
92 G4double PHEP2; // py in GeV
93 G4double PHEP3; // pz in GeV
94 G4double PHEP5; // mass in GeV
95
96 inputFile >> ISTHEP >> IDHEP >> JDAHEP1 >> JDAHEP2
97 >> PHEP1 >> PHEP2 >> PHEP3 >> PHEP5;
98
99 // create G4PrimaryParticle object
100 G4PrimaryParticle* particle
101 = new G4PrimaryParticle( IDHEP, PHEP1*GeV, PHEP2*GeV, PHEP3*GeV );
102 particle->SetMass( PHEP5*GeV );
103
104 // create G4HEPEvtParticle object
105 G4HEPEvtParticle* hepParticle
106 = new G4HEPEvtParticle( particle, ISTHEP, JDAHEP1, JDAHEP2 );
107
108 // Store
109 HPlist.push_back( hepParticle );
110 }
111
112 // check if there is at least one particle
113 if( HPlist.size() == 0 ) return;
114
115 // make connection between daughter particles decayed from
116 // the same mother
117 for( size_t i=0; i<HPlist.size(); i++ )
118 {
119 if( HPlist[i]->GetJDAHEP1() > 0 ) // it has daughters
120 {
121 G4int jda1 = HPlist[i]->GetJDAHEP1()-1; // FORTRAN index starts from 1
122 G4int jda2 = HPlist[i]->GetJDAHEP2()-1; // but C++ starts from 0.
123 G4PrimaryParticle* mother = HPlist[i]->GetTheParticle();
124 for( G4int j=jda1; j<=jda2; j++ )
125 {
126 G4PrimaryParticle* daughter = HPlist[j]->GetTheParticle();
127 if(HPlist[j]->GetISTHEP()>0)
128 {
129 mother->SetDaughter( daughter );
130 HPlist[j]->Done();
131 }
132 }
133 }
134 }
135
136 // create G4PrimaryVertex object
137 G4PrimaryVertex* vertex = new G4PrimaryVertex(particle_position,particle_time);
138
139 // put initial particles to the vertex
140 for( size_t ii=0; ii<HPlist.size(); ii++ )
141 {
142 if( HPlist[ii]->GetISTHEP() > 0 ) // ISTHEP of daughters had been
143 // set to negative
144 {
145 G4PrimaryParticle* initialParticle = HPlist[ii]->GetTheParticle();
146 vertex->SetPrimary( initialParticle );
147 }
148 }
149
150 // clear G4HEPEvtParticles
151 //HPlist.clearAndDestroy();
152 for(size_t iii=0;iii<HPlist.size();iii++)
153 { delete HPlist[iii]; }
154 HPlist.clear();
155
156 // Put the vertex to G4Event object
157 evt->AddPrimaryVertex( vertex );
158}
159
Note: See TracBrowser for help on using the repository browser.