G4HEPEvtInterface.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 //
00027 // $Id$
00028 //
00029 // 
00030 // --------------------------------------------------------------------
00031 
00032 #include "G4HEPEvtInterface.hh"
00033 
00034 #include "G4Types.hh"
00035 #include "G4SystemOfUnits.hh"
00036 
00037 #include "G4ios.hh"
00038 #include "G4PrimaryVertex.hh"
00039 #include "G4PrimaryParticle.hh"
00040 #include "G4HEPEvtParticle.hh"
00041 #include "G4Event.hh"
00042 
00043 G4HEPEvtInterface::G4HEPEvtInterface(char* evfile)
00044 {
00045   inputFile.open(evfile);
00046   if (inputFile) {
00047     fileName = evfile;
00048   }
00049   else {
00050     G4Exception("G4HEPEvtInterface::G4HEPEvtInterface","Event0201",FatalException,
00051     "G4HEPEvtInterface:: cannot open file.");
00052   }
00053   G4ThreeVector zero;
00054   particle_position = zero;
00055   particle_time = 0.0;
00056 
00057 }
00058 
00059 G4HEPEvtInterface::G4HEPEvtInterface(G4String evfile)
00060 {
00061   const char* fn = evfile.data();
00062   inputFile.open((char*)fn);
00063   if (inputFile) {
00064     fileName = evfile;
00065   }
00066   else {
00067     G4Exception("G4HEPEvtInterface::G4HEPEvtInterface","Event0201",FatalException,
00068     "G4HEPEvtInterface:: cannot open file.");
00069   }
00070   G4ThreeVector zero;
00071   particle_position = zero;
00072   particle_time = 0.0;
00073 }
00074 
00075 G4HEPEvtInterface::~G4HEPEvtInterface()
00076 {;}
00077 
00078 void G4HEPEvtInterface::GeneratePrimaryVertex(G4Event* evt)
00079 {
00080   G4int NHEP;  // number of entries
00081   inputFile >> NHEP;
00082   if( inputFile.eof() ) 
00083   {
00084     G4Exception("G4HEPEvtInterface::GeneratePrimaryVertex","Event0202",
00085     JustWarning,"End-Of-File : HEPEvt input file");
00086     return;
00087   }
00088 
00089   for( G4int IHEP=0; IHEP<NHEP; IHEP++ )
00090   {
00091     G4int ISTHEP;   // status code
00092     G4int IDHEP;    // PDG code
00093     G4int JDAHEP1;  // first daughter
00094     G4int JDAHEP2;  // last daughter
00095     G4double PHEP1; // px in GeV
00096     G4double PHEP2; // py in GeV
00097     G4double PHEP3; // pz in GeV
00098     G4double PHEP5; // mass in GeV
00099 
00100     inputFile >> ISTHEP >> IDHEP >> JDAHEP1 >> JDAHEP2
00101        >> PHEP1 >> PHEP2 >> PHEP3 >> PHEP5;
00102 
00103     // create G4PrimaryParticle object
00104     G4PrimaryParticle* particle 
00105       = new G4PrimaryParticle( IDHEP );
00106     particle->SetMass( PHEP5*GeV );
00107     particle->SetMomentum(PHEP1*GeV, PHEP2*GeV, PHEP3*GeV );
00108 
00109     // create G4HEPEvtParticle object
00110     G4HEPEvtParticle* hepParticle
00111       = new G4HEPEvtParticle( particle, ISTHEP, JDAHEP1, JDAHEP2 );
00112 
00113     // Store
00114     HPlist.push_back( hepParticle );
00115   }
00116 
00117   // check if there is at least one particle
00118   if( HPlist.size() == 0 ) return; 
00119 
00120   // make connection between daughter particles decayed from 
00121   // the same mother
00122   for( size_t i=0; i<HPlist.size(); i++ )
00123   {
00124     if( HPlist[i]->GetJDAHEP1() > 0 ) //  it has daughters
00125     {
00126       G4int jda1 = HPlist[i]->GetJDAHEP1()-1; // FORTRAN index starts from 1
00127       G4int jda2 = HPlist[i]->GetJDAHEP2()-1; // but C++ starts from 0.
00128       G4PrimaryParticle* mother = HPlist[i]->GetTheParticle();
00129       for( G4int j=jda1; j<=jda2; j++ )
00130       {
00131         G4PrimaryParticle* daughter = HPlist[j]->GetTheParticle();
00132         if(HPlist[j]->GetISTHEP()>0)
00133         {
00134           mother->SetDaughter( daughter );
00135           HPlist[j]->Done();
00136         }
00137       }
00138     }
00139   }
00140 
00141   // create G4PrimaryVertex object
00142   G4PrimaryVertex* vertex = new G4PrimaryVertex(particle_position,particle_time);
00143 
00144   // put initial particles to the vertex
00145   for( size_t ii=0; ii<HPlist.size(); ii++ )
00146   {
00147     if( HPlist[ii]->GetISTHEP() > 0 ) // ISTHEP of daughters had been 
00148                                        // set to negative
00149     {
00150       G4PrimaryParticle* initialParticle = HPlist[ii]->GetTheParticle();
00151       vertex->SetPrimary( initialParticle );
00152     }
00153   }
00154 
00155   // clear G4HEPEvtParticles
00156   //HPlist.clearAndDestroy();
00157   for(size_t iii=0;iii<HPlist.size();iii++)
00158   { delete HPlist[iii]; }
00159   HPlist.clear();
00160 
00161   // Put the vertex to G4Event object
00162   evt->AddPrimaryVertex( vertex );
00163 }
00164 

Generated on Mon May 27 17:48:29 2013 for Geant4 by  doxygen 1.4.7