Geant4 Cross Reference |
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 /// \file eventgenerator/HepMC/HepMCEx01/src/HepMCG4Interface.cc 27 /// \brief Implementation of the HepMCG4Interface class 28 // 29 // 30 31 #include "HepMCG4Interface.hh" 32 33 #include "G4Event.hh" 34 #include "G4LorentzVector.hh" 35 #include "G4PhysicalConstants.hh" 36 #include "G4PrimaryParticle.hh" 37 #include "G4PrimaryVertex.hh" 38 #include "G4RunManager.hh" 39 #include "G4SystemOfUnits.hh" 40 #include "G4TransportationManager.hh" 41 42 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 43 HepMCG4Interface::HepMCG4Interface() : hepmcEvent(0) {} 44 45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 46 HepMCG4Interface::~HepMCG4Interface() 47 { 48 delete hepmcEvent; 49 } 50 51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 52 G4bool HepMCG4Interface::CheckVertexInsideWorld(const G4ThreeVector& pos) const 53 { 54 G4Navigator* navigator = 55 G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking(); 56 57 G4VPhysicalVolume* world = navigator->GetWorldVolume(); 58 G4VSolid* solid = world->GetLogicalVolume()->GetSolid(); 59 EInside qinside = solid->Inside(pos); 60 61 if (qinside != kInside) 62 return false; 63 else 64 return true; 65 } 66 67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 68 void HepMCG4Interface::HepMC2G4(const HepMC::GenEvent* hepmcevt, G4Event* g4event) 69 { 70 for (HepMC::GenEvent::vertex_const_iterator vitr = hepmcevt->vertices_begin(); 71 vitr != hepmcevt->vertices_end(); ++vitr) 72 { // loop for vertex ... 73 74 // real vertex? 75 G4bool qvtx = false; 76 for (HepMC::GenVertex::particle_iterator pitr = (*vitr)->particles_begin(HepMC::children); 77 pitr != (*vitr)->particles_end(HepMC::children); ++pitr) 78 { 79 if (!(*pitr)->end_vertex() && (*pitr)->status() == 1) { 80 qvtx = true; 81 break; 82 } 83 } 84 if (!qvtx) continue; 85 86 // check world boundary 87 HepMC::FourVector pos = (*vitr)->position(); 88 G4LorentzVector xvtx(pos.x(), pos.y(), pos.z(), pos.t()); 89 if (!CheckVertexInsideWorld(xvtx.vect() * mm)) continue; 90 91 // create G4PrimaryVertex and associated G4PrimaryParticles 92 G4PrimaryVertex* g4vtx = 93 new G4PrimaryVertex(xvtx.x() * mm, xvtx.y() * mm, xvtx.z() * mm, xvtx.t() * mm / c_light); 94 95 for (HepMC::GenVertex::particle_iterator vpitr = (*vitr)->particles_begin(HepMC::children); 96 vpitr != (*vitr)->particles_end(HepMC::children); ++vpitr) 97 { 98 if ((*vpitr)->status() != 1) continue; 99 100 G4int pdgcode = (*vpitr)->pdg_id(); 101 pos = (*vpitr)->momentum(); 102 G4LorentzVector p(pos.px(), pos.py(), pos.pz(), pos.e()); 103 G4PrimaryParticle* g4prim = 104 new G4PrimaryParticle(pdgcode, p.x() * GeV, p.y() * GeV, p.z() * GeV); 105 106 g4vtx->SetPrimary(g4prim); 107 } 108 g4event->AddPrimaryVertex(g4vtx); 109 } 110 } 111 112 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 113 HepMC::GenEvent* HepMCG4Interface::GenerateHepMCEvent() 114 { 115 HepMC::GenEvent* aevent = new HepMC::GenEvent(); 116 return aevent; 117 } 118 119 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo...... 120 void HepMCG4Interface::GeneratePrimaryVertex(G4Event* anEvent) 121 { 122 // delete previous event object 123 delete hepmcEvent; 124 125 // generate next event 126 hepmcEvent = GenerateHepMCEvent(); 127 if (!hepmcEvent) { 128 G4cout << "HepMCInterface: no generated particles. run terminated..." << G4endl; 129 G4RunManager::GetRunManager()->AbortRun(); 130 return; 131 } 132 HepMC2G4(hepmcEvent, anEvent); 133 } 134