File indexing completed on 2026-05-13 08:05:19
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026
0027
0028
0029 #include "HepMCG4Interface.hh"
0030
0031 #include "G4Event.hh"
0032 #include "G4LorentzVector.hh"
0033 #include "G4PhysicalConstants.hh"
0034 #include "G4PrimaryParticle.hh"
0035 #include "G4PrimaryVertex.hh"
0036 #include "G4RunManager.hh"
0037 #include "G4SystemOfUnits.hh"
0038 #include "G4TransportationManager.hh"
0039
0040
0041 HepMCG4Interface::HepMCG4Interface() : hepmcEvent(0) {}
0042
0043
0044 HepMCG4Interface::~HepMCG4Interface()
0045 {
0046 delete hepmcEvent;
0047 }
0048
0049
0050 G4bool HepMCG4Interface::CheckVertexInsideWorld(const G4ThreeVector& pos) const
0051 {
0052 G4Navigator* navigator =
0053 G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking();
0054
0055 G4VPhysicalVolume* world = navigator->GetWorldVolume();
0056 G4VSolid* solid = world->GetLogicalVolume()->GetSolid();
0057 EInside qinside = solid->Inside(pos);
0058
0059 if (qinside != kInside)
0060 return false;
0061 else
0062 return true;
0063 }
0064
0065
0066 void HepMCG4Interface::HepMC2G4(const HepMC::GenEvent* hepmcevt, G4Event* g4event)
0067 {
0068 for (HepMC::GenEvent::vertex_const_iterator vitr = hepmcevt->vertices_begin();
0069 vitr != hepmcevt->vertices_end(); ++vitr)
0070 {
0071
0072
0073 G4bool qvtx = false;
0074 for (HepMC::GenVertex::particle_iterator pitr = (*vitr)->particles_begin(HepMC::children);
0075 pitr != (*vitr)->particles_end(HepMC::children); ++pitr)
0076 {
0077 if (!(*pitr)->end_vertex() && (*pitr)->status() == 1) {
0078 qvtx = true;
0079 break;
0080 }
0081 }
0082 if (!qvtx) continue;
0083
0084
0085 HepMC::FourVector pos = (*vitr)->position();
0086 G4LorentzVector xvtx(pos.x(), pos.y(), pos.z(), pos.t());
0087 if (!CheckVertexInsideWorld(xvtx.vect() * mm)) continue;
0088
0089
0090 G4PrimaryVertex* g4vtx =
0091 new G4PrimaryVertex(xvtx.x() * mm, xvtx.y() * mm, xvtx.z() * mm, xvtx.t() * mm / c_light);
0092
0093 for (HepMC::GenVertex::particle_iterator vpitr = (*vitr)->particles_begin(HepMC::children);
0094 vpitr != (*vitr)->particles_end(HepMC::children); ++vpitr)
0095 {
0096 if ((*vpitr)->status() != 1) continue;
0097
0098 G4int pdgcode = (*vpitr)->pdg_id();
0099 pos = (*vpitr)->momentum();
0100 G4LorentzVector p(pos.px(), pos.py(), pos.pz(), pos.e());
0101 G4PrimaryParticle* g4prim =
0102 new G4PrimaryParticle(pdgcode, p.x() * GeV, p.y() * GeV, p.z() * GeV);
0103
0104 g4vtx->SetPrimary(g4prim);
0105 }
0106 g4event->AddPrimaryVertex(g4vtx);
0107 }
0108 }
0109
0110
0111 HepMC::GenEvent* HepMCG4Interface::GenerateHepMCEvent()
0112 {
0113 HepMC::GenEvent* aevent = new HepMC::GenEvent();
0114 return aevent;
0115 }
0116
0117
0118 void HepMCG4Interface::GeneratePrimaryVertex(G4Event* anEvent)
0119 {
0120
0121 delete hepmcEvent;
0122
0123
0124 hepmcEvent = GenerateHepMCEvent();
0125 if (!hepmcEvent) {
0126 G4cout << "HepMCInterface: no generated particles. run terminated..." << G4endl;
0127 G4RunManager::GetRunManager()->AbortRun();
0128 return;
0129 }
0130 HepMC2G4(hepmcEvent, anEvent);
0131 }