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